Changeset 3983 for trunk/LMDZ.COMMON/libf/evolution/read_XIOS_data.F90
- Timestamp:
- Dec 8, 2025, 11:27:43 AM (7 days ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.COMMON/libf/evolution/read_XIOS_data.F90
r3977 r3983 19 19 contains 20 20 !======================================================================= 21 SUBROUTINE read_PCM_data(ngrid,nslope,nsoil_PCM,nsol,ps_avg,tsurf_avg,tsurf_avg_y1,tsoil_avg,tsoil_ts,watersurf_density_avg,d_h2oice,d_co2ice, & 22 ps_ts,q_h2o_ts,q_co2_ts,watersoil_density_ts) 23 24 use compute_tend_mod, only: compute_tend 21 SUBROUTINE read_PCM_data(ngrid,nslope,nsoil_PCM,nsol,h2ofrost_PCM,co2frost_PCM,ps_avg,tsurf_avg,tsurf_avg_y1,tsoil_avg,tsoil_ts,watersurf_density_avg,d_h2oice,d_co2ice, & 22 ps_ts,q_h2o_ts,q_co2_ts,watersoil_density_ts,min_h2oice,min_co2ice) 23 25 24 use grid_conversion, only: lonlat2vect 26 25 use comsoil_h_PEM, only: soil_pem 26 use compute_tend_mod, only: compute_tend 27 use metamorphism, only: compute_frost 27 28 28 29 implicit none … … 30 31 ! Arguments 31 32 !---------- 32 integer, intent(in) :: ngrid, nslope, nsoil_PCM, nsol 33 integer, intent(in) :: ngrid, nslope, nsoil_PCM, nsol 34 real, dimension(ngrid,nslope), intent(in) :: h2ofrost_PCM, co2frost_PCM 33 35 real, dimension(ngrid), intent(out) :: ps_avg 34 real, dimension(ngrid,nslope), intent(out) :: tsurf_avg, tsurf_avg_y1, watersurf_density_avg, d_h2oice, d_co2ice 36 real, dimension(ngrid,nslope), intent(out) :: tsurf_avg, tsurf_avg_y1, watersurf_density_avg, d_h2oice, d_co2ice, min_h2oice, min_co2ice 35 37 real, dimension(ngrid,nsoil_PCM,nslope), intent(out) :: tsoil_avg 36 38 real, dimension(ngrid,nsol), intent(out) :: ps_ts, q_h2o_ts, q_co2_ts … … 44 46 real, dimension(:,:,:,:), allocatable :: var_read_4d 45 47 character(:), allocatable :: num ! For reading slope variables 46 real, dimension(ngrid,nslope,2) :: min_h2o ice, min_co2ice, min_h2ofrost, min_co2frost48 real, dimension(ngrid,nslope,2) :: min_h2operice, min_co2perice, min_h2ofrost, min_co2frost 47 49 48 50 ! Code 49 51 !----- 50 52 ! Initialization 51 min_h2o ice = 0.52 min_co2 ice = 0.53 min_h2operice = 0. 54 min_co2perice = 0. 53 55 min_h2ofrost = 0. 54 56 min_co2frost = 0. … … 80 82 call get_var('co2ice'//num,var_read_2d) ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,min_co2frost(:,islope,1)) 81 83 call get_var('h2o_ice_s'//num,var_read_2d) ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,min_h2ofrost(:,islope,1)) 82 call get_var('watercap'//num,var_read_2d) ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,min_h2o ice(:,islope,1))83 call get_var('perennial_co2ice'//num,var_read_2d); call lonlat2vect(nlon,nlat,ngrid,var_read_2d,min_co2 ice(:,islope,1))84 call get_var('watercap'//num,var_read_2d) ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,min_h2operice(:,islope,1)) 85 call get_var('perennial_co2ice'//num,var_read_2d); call lonlat2vect(nlon,nlat,ngrid,var_read_2d,min_co2perice(:,islope,1)) 84 86 call get_var('tsurf'//num,var_read_2d) ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,tsurf_avg_y1(:,islope)) 85 87 enddo … … 95 97 96 98 ! Allocate and read the variables 97 call get_var('ps',var_read_2d) ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,ps_avg)99 call get_var('ps',var_read_2d); call lonlat2vect(nlon,nlat,ngrid,var_read_2d,ps_avg) 98 100 do islope = 1,nslope 99 101 if (nslope /= 1) then … … 105 107 call get_var('co2ice'//num,var_read_2d) ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,min_co2frost(:,islope,2)) 106 108 call get_var('h2o_ice_s'//num,var_read_2d) ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,min_h2ofrost(:,islope,2)) 107 call get_var('watercap'//num,var_read_2d) ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,min_h2o ice(:,islope,2))108 call get_var('perennial_co2ice'//num,var_read_2d); call lonlat2vect(nlon,nlat,ngrid,var_read_2d,min_co2 ice(:,islope,2))109 call get_var('watercap'//num,var_read_2d) ; call lonlat2vect(nlon,nlat,ngrid,var_read_2d,min_h2operice(:,islope,2)) 110 call get_var('perennial_co2ice'//num,var_read_2d); call lonlat2vect(nlon,nlat,ngrid,var_read_2d,min_co2perice(:,islope,2)) 109 111 if (soil_pem) then 110 112 call get_var('soiltemp'//num,var_read_3d) … … 167 169 168 170 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 171 ! Compute frost from yearly minima 172 call compute_frost(ngrid,nslope,h2ofrost_PCM,min_h2ofrost,co2frost_PCM,min_co2frost) 173 174 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 169 175 ! Compute ice tendencies from yearly minima 170 176 write(*,*) '> Computing surface ice tendencies' 171 call compute_tend(ngrid,nslope,min_h2o frost + min_h2oice,d_h2oice)177 call compute_tend(ngrid,nslope,min_h2operice + min_h2ofrost,d_h2oice) 172 178 write(*,*) 'H2O ice tendencies (min/max):', minval(d_h2oice), maxval(d_h2oice) 173 call compute_tend(ngrid,nslope,min_co2 frost + min_co2ice,d_co2ice)179 call compute_tend(ngrid,nslope,min_co2perice + min_co2frost,d_co2ice) 174 180 write(*,*) 'CO2 ice tendencies (min/max):', minval(d_co2ice), maxval(d_co2ice) 175 181
Note: See TracChangeset
for help on using the changeset viewer.
