| 1 | MODULE read_data_PCM_mod |
|---|
| 2 | |
|---|
| 3 | use netcdf, only: nf90_open, NF90_NOWRITE, nf90_noerr, nf90_strerror, & |
|---|
| 4 | nf90_get_var, nf90_inq_varid, nf90_inq_dimid, & |
|---|
| 5 | nf90_inquire_dimension, nf90_close |
|---|
| 6 | |
|---|
| 7 | implicit none |
|---|
| 8 | |
|---|
| 9 | character(13), parameter :: modname = 'read_data_PCM' |
|---|
| 10 | character(256) :: msg ! for reading |
|---|
| 11 | integer :: fID, vID ! for reading |
|---|
| 12 | |
|---|
| 13 | !======================================================================= |
|---|
| 14 | contains |
|---|
| 15 | !======================================================================= |
|---|
| 16 | |
|---|
| 17 | SUBROUTINE read_data_PCM(filename1,filename2,timelen,iim_input,jjm_input,ngrid,nslope,vmr_co2_PCM,ps_timeseries,ps_avg,tsurf_avg_yr1,tsurf_avg, & |
|---|
| 18 | tsoil_avg,tsoil_timeseries,min_co2_ice,min_h2o_ice,q_co2,q_h2o,watersurf_density_avg,watersoil_density_timeseries) |
|---|
| 19 | use comsoil_h, only: nsoilmx |
|---|
| 20 | use comsoil_h_PEM, only: soil_pem |
|---|
| 21 | use constants_marspem_mod, only: m_co2, m_noco2 |
|---|
| 22 | |
|---|
| 23 | implicit none |
|---|
| 24 | |
|---|
| 25 | !======================================================================= |
|---|
| 26 | ! |
|---|
| 27 | ! Purpose: Read initial confitions file from the PCM |
|---|
| 28 | ! |
|---|
| 29 | ! Authors: JBC |
|---|
| 30 | !======================================================================= |
|---|
| 31 | |
|---|
| 32 | include "dimensions.h" |
|---|
| 33 | |
|---|
| 34 | !======================================================================= |
|---|
| 35 | ! Inputs: |
|---|
| 36 | character(*), intent(in) :: filename1, filename2 ! File names |
|---|
| 37 | integer, intent(in) :: timelen ! Number of times stored in the file |
|---|
| 38 | integer, intent(in) :: iim_input, jjm_input, ngrid, nslope ! Number of points in the lat x lon dynamical grid, number of subgrid slopes |
|---|
| 39 | ! Ouputs |
|---|
| 40 | real, dimension(ngrid,nslope,2), intent(out) :: min_co2_ice ! Minimum of co2 ice per slope of the year [kg/m^2] |
|---|
| 41 | real, dimension(ngrid,nslope,2), intent(out) :: min_h2o_ice ! Minimum of h2o ice per slope of the year [kg/m^2] |
|---|
| 42 | real, dimension(ngrid,timelen), intent(out) :: vmr_co2_PCM ! Grid points x Times co2 volume mixing ratio retrieve from the PCM [m^3/m^3] |
|---|
| 43 | real, dimension(ngrid,timelen), intent(out) :: q_co2 ! CO2 mass mixing ratio in the first layer [kg/m^3] |
|---|
| 44 | real, dimension(ngrid,timelen), intent(out) :: q_h2o ! H2O mass mixing ratio in the first layer [kg/m^3] |
|---|
| 45 | real, dimension(ngrid,timelen), intent(out) :: ps_timeseries ! Surface pressure timeseries [Pa] |
|---|
| 46 | real, dimension(ngrid), intent(out) :: ps_avg ! Averaged surface pressure [K] |
|---|
| 47 | real, dimension(ngrid,nslope), intent(out) :: tsurf_avg ! Averaged surface temperature [K] |
|---|
| 48 | real, dimension(ngrid,nslope), intent(out) :: tsurf_avg_yr1 ! Averaged surface temperature for year 1 [K] |
|---|
| 49 | real, dimension(ngrid,nsoilmx,nslope), intent(out) :: tsoil_avg ! Averaged soil temperature for year 2 [K] |
|---|
| 50 | real, dimension(ngrid,nsoilmx,nslope,timelen), intent(out) :: tsoil_timeseries ! Soil temperature timeseries [K] |
|---|
| 51 | real, dimension(ngrid,nslope), intent(out) :: watersurf_density_avg ! Water density at the surface [kg/m^3] |
|---|
| 52 | real, dimension(ngrid,nsoilmx,nslope,timelen), intent(out) :: watersoil_density_timeseries ! Water density timeseries in the soil layer [kg/m^3] |
|---|
| 53 | ! Local variables |
|---|
| 54 | integer :: i, j, l, islope ! Loop variables |
|---|
| 55 | real :: A, B ! Intermediate variables to compute the mean molar mass of the layer |
|---|
| 56 | character(2) :: num ! For reading sloped variables |
|---|
| 57 | |
|---|
| 58 | real, dimension(:,:), allocatable :: var2_read ! Variables for reading (2 dimensions) |
|---|
| 59 | real, dimension(:,:,:), allocatable :: var3_read_1, var3_read_2, var3_read_3 ! Variables for reading (3 dimensions) |
|---|
| 60 | real, dimension(:,:,:,:), allocatable :: var4_read ! Variables for reading (4 dimensions) |
|---|
| 61 | !----------------------------------------------------------------------- |
|---|
| 62 | |
|---|
| 63 | !------------------------------- Year 1 -------------------------------- |
|---|
| 64 | ! Open the NetCDF file |
|---|
| 65 | write(*,*) "Opening "//filename1//"..." |
|---|
| 66 | call error_msg(NF90_OPEN(filename1,NF90_NOWRITE,fID),"open",filename1) |
|---|
| 67 | |
|---|
| 68 | ! Download the data from the file |
|---|
| 69 | allocate(var2_read(iim_input + 1,jjm_input + 1),var3_read_1(iim_input + 1,jjm_input + 1,timelen),var3_read_2(iim_input + 1,jjm_input + 1,timelen)) |
|---|
| 70 | |
|---|
| 71 | if (nslope == 1) then ! There is no slope |
|---|
| 72 | ! CO2 ice |
|---|
| 73 | !-------- |
|---|
| 74 | call get_var3("co2ice",var3_read_1) |
|---|
| 75 | where (var3_read_1 < 0.) var3_read_1 = 0. |
|---|
| 76 | write(*,*) "Data for co2_ice downloaded." |
|---|
| 77 | #ifndef CPP_STD |
|---|
| 78 | call get_var3("perennial_co2ice",var3_read_2) |
|---|
| 79 | write(*,*) "Data for perennial_co2ice downloaded." |
|---|
| 80 | #endif |
|---|
| 81 | |
|---|
| 82 | ! Compute the minimum over the year for each point |
|---|
| 83 | var2_read = minval(var3_read_1 + var3_read_2,3) |
|---|
| 84 | #ifndef CPP_1D |
|---|
| 85 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_co2_ice(:,1,1)) |
|---|
| 86 | #else |
|---|
| 87 | min_co2_ice(1,1,1) = var2_read(1,1) |
|---|
| 88 | #endif |
|---|
| 89 | write(*,*) "Min of co2_ice computed." |
|---|
| 90 | |
|---|
| 91 | ! H2O ice |
|---|
| 92 | !-------- |
|---|
| 93 | call get_var3("h2o_ice_s",var3_read_1) |
|---|
| 94 | where (var3_read_1 < 0.) var3_read_1 = 0. |
|---|
| 95 | write(*,*) "Data for h2o_ice_s downloaded." |
|---|
| 96 | #ifndef CPP_STD |
|---|
| 97 | call get_var3("watercap",var3_read_2) |
|---|
| 98 | write(*,*) "Data for watercap downloaded." |
|---|
| 99 | #endif |
|---|
| 100 | |
|---|
| 101 | ! Compute the minimum over the year for each point |
|---|
| 102 | var2_read = minval(var3_read_1 + var3_read_2,3) |
|---|
| 103 | #ifndef CPP_1D |
|---|
| 104 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_h2o_ice(:,1,1)) |
|---|
| 105 | #else |
|---|
| 106 | min_h2o_ice(1,1,1) = var2_read(1,1) |
|---|
| 107 | #endif |
|---|
| 108 | write(*,*) "Min of h2o_ice computed." |
|---|
| 109 | |
|---|
| 110 | ! Tsurf |
|---|
| 111 | !------ |
|---|
| 112 | call get_var3("tsurf",var3_read_1) |
|---|
| 113 | write(*,*) "Data for tsurf downloaded." |
|---|
| 114 | |
|---|
| 115 | ! Compute average over the year for each point |
|---|
| 116 | var2_read = sum(var3_read_1,3)/timelen |
|---|
| 117 | #ifndef CPP_1D |
|---|
| 118 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,tsurf_avg_yr1(:,1)) |
|---|
| 119 | #else |
|---|
| 120 | tsurf_avg_yr1(1,1) = var2_read(1,1) |
|---|
| 121 | #endif |
|---|
| 122 | write(*,*) "Average of tsurf computed." |
|---|
| 123 | else ! We use slopes |
|---|
| 124 | do islope = 1,nslope |
|---|
| 125 | write(num,'(i2.2)') islope |
|---|
| 126 | |
|---|
| 127 | ! CO2 ice |
|---|
| 128 | !-------- |
|---|
| 129 | call get_var3("co2ice_slope"//num,var3_read_1) |
|---|
| 130 | where (var3_read_1 < 0.) var3_read_1 = 0. |
|---|
| 131 | write(*,*) "Data for co2_ice_slope"//num//" downloaded." |
|---|
| 132 | #ifndef CPP_STD |
|---|
| 133 | call get_var3("perennial_co2ice_slope"//num,var3_read_2) |
|---|
| 134 | write(*,*) "Data for perennial_co2ice_slope"//num//" downloaded." |
|---|
| 135 | #endif |
|---|
| 136 | |
|---|
| 137 | ! Compute the minimum over the year for each point |
|---|
| 138 | var2_read = minval(var3_read_1 + var3_read_2,3) |
|---|
| 139 | #ifndef CPP_1D |
|---|
| 140 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_co2_ice(:,islope,1)) |
|---|
| 141 | #else |
|---|
| 142 | min_co2_ice(1,islope,1) = var2_read(1,1) |
|---|
| 143 | #endif |
|---|
| 144 | write(*,*) "Min of co2_ice computed for slope "//num//"." |
|---|
| 145 | |
|---|
| 146 | ! H2O ice |
|---|
| 147 | !-------- |
|---|
| 148 | call get_var3("h2o_ice_s_slope"//num,var3_read_1) |
|---|
| 149 | where (var3_read_1 < 0.) var3_read_1 = 0. |
|---|
| 150 | write(*,*) "Data for h2o_ice_s_slope"//num//" downloaded." |
|---|
| 151 | #ifndef CPP_STD |
|---|
| 152 | call get_var3("watercap_slope"//num,var3_read_2) |
|---|
| 153 | write(*,*) "Data for watercap_slope"//num//" downloaded." |
|---|
| 154 | #endif |
|---|
| 155 | |
|---|
| 156 | ! Compute the minimum over the year for each point |
|---|
| 157 | var2_read = minval(var3_read_1 + var3_read_2,3) |
|---|
| 158 | #ifndef CPP_1D |
|---|
| 159 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_h2o_ice(:,islope,1)) |
|---|
| 160 | #else |
|---|
| 161 | min_h2o_ice(1,islope,1) = var2_read(1,1) |
|---|
| 162 | #endif |
|---|
| 163 | write(*,*) "Min of h2o_ice computed for slope "//num//"." |
|---|
| 164 | |
|---|
| 165 | ! Tsurf |
|---|
| 166 | !------ |
|---|
| 167 | call get_var3("tsurf_slope"//num,var3_read_1) |
|---|
| 168 | write(*,*) "Data for tsurf_slope"//num//" downloaded." |
|---|
| 169 | |
|---|
| 170 | ! Compute average over the year for each point |
|---|
| 171 | var2_read = sum(var3_read_1,3)/timelen |
|---|
| 172 | #ifndef CPP_1D |
|---|
| 173 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,tsurf_avg_yr1(:,islope)) |
|---|
| 174 | #else |
|---|
| 175 | tsurf_avg_yr1(1,islope) = var2_read(1,1) |
|---|
| 176 | #endif |
|---|
| 177 | write(*,*) "Average of tsurf computed for slope "//num//"." |
|---|
| 178 | enddo |
|---|
| 179 | endif |
|---|
| 180 | |
|---|
| 181 | ! Close the NetCDF file |
|---|
| 182 | call error_msg(nf90_close(fID),"close",filename1) |
|---|
| 183 | write(*,*) '> "'//filename1//'" processed!' |
|---|
| 184 | |
|---|
| 185 | !------------------------------- Year 2 -------------------------------- |
|---|
| 186 | ! Open the NetCDF file |
|---|
| 187 | write(*,*) "Opening "//filename2//"..." |
|---|
| 188 | call error_msg(NF90_OPEN(filename2,NF90_NOWRITE,fID),"open",filename2) |
|---|
| 189 | |
|---|
| 190 | ! Download the data from the file |
|---|
| 191 | allocate(var3_read_3(iim_input + 1,jjm_input + 1,nsoilmx),var4_read(iim_input + 1,jjm_input + 1,nsoilmx,timelen)) |
|---|
| 192 | |
|---|
| 193 | ! Surface pressure |
|---|
| 194 | !----------------- |
|---|
| 195 | call get_var3("ps",var3_read_1) |
|---|
| 196 | write(*,*) "Data for surface pressure downloaded." |
|---|
| 197 | |
|---|
| 198 | ! Compute average over the year for each point |
|---|
| 199 | var2_read = sum(var3_read_1,3)/timelen |
|---|
| 200 | #ifndef CPP_1D |
|---|
| 201 | call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,var3_read_1,ps_timeseries) |
|---|
| 202 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,ps_avg) |
|---|
| 203 | #else |
|---|
| 204 | ps_timeseries(1,:) = var3_read_1(1,:) |
|---|
| 205 | ps_avg(1) = var2_read(1,1) |
|---|
| 206 | #endif |
|---|
| 207 | write(*,*) "Average of surface pressure computed." |
|---|
| 208 | |
|---|
| 209 | ! CO2 vmr |
|---|
| 210 | !-------- |
|---|
| 211 | call get_var3("co2_layer1",var3_read_1) |
|---|
| 212 | where (var3_read_1 < 0.) |
|---|
| 213 | var3_read_1 = 1.e-10 |
|---|
| 214 | else where (var3_read_1 > 1.) |
|---|
| 215 | var3_read_1 = 1. |
|---|
| 216 | end where |
|---|
| 217 | #ifndef CPP_1D |
|---|
| 218 | call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,var3_read_1,q_co2) |
|---|
| 219 | #else |
|---|
| 220 | q_co2(1,:) = var3_read_1(1,1,:) |
|---|
| 221 | #endif |
|---|
| 222 | A = (1./m_co2 - 1./m_noco2) |
|---|
| 223 | B = 1./m_noco2 |
|---|
| 224 | vmr_co2_PCM = q_co2/(A*q_co2 + B)/m_co2 |
|---|
| 225 | write(*,*) "Data for CO2 vmr downloaded." |
|---|
| 226 | |
|---|
| 227 | ! H2O vmr |
|---|
| 228 | !-------- |
|---|
| 229 | call get_var3("h2o_layer1",var3_read_1) |
|---|
| 230 | where (var3_read_1 < 0.) |
|---|
| 231 | var3_read_1 = 1.e-10 |
|---|
| 232 | else where (var3_read_1 > 1.) |
|---|
| 233 | var3_read_1 = 1. |
|---|
| 234 | end where |
|---|
| 235 | #ifndef CPP_1D |
|---|
| 236 | call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,var3_read_1,q_h2o) |
|---|
| 237 | #else |
|---|
| 238 | q_h2o(1,:) = var3_read_1(1,1,:) |
|---|
| 239 | #endif |
|---|
| 240 | write(*,*) "Data for H2O vmr downloaded." |
|---|
| 241 | |
|---|
| 242 | if (nslope == 1) then ! There is no slope |
|---|
| 243 | ! CO2 ice |
|---|
| 244 | !-------- |
|---|
| 245 | call get_var3("co2ice",var3_read_1) |
|---|
| 246 | where (var3_read_1 < 0.) var3_read_1 = 0. |
|---|
| 247 | write(*,*) "Data for co2_ice downloaded." |
|---|
| 248 | #ifndef CPP_STD |
|---|
| 249 | call get_var3("perennial_co2ice",var3_read_2) |
|---|
| 250 | write(*,*) "Data for perennial_co2ice downloaded." |
|---|
| 251 | #endif |
|---|
| 252 | |
|---|
| 253 | ! Compute the minimum over the year for each point |
|---|
| 254 | var2_read = minval(var3_read_1 + var3_read_2,3) |
|---|
| 255 | #ifndef CPP_1D |
|---|
| 256 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_co2_ice(:,1,1)) |
|---|
| 257 | #else |
|---|
| 258 | min_co2_ice(1,1,1) = var2_read(1,1) |
|---|
| 259 | #endif |
|---|
| 260 | write(*,*) "Min of co2_ice computed." |
|---|
| 261 | |
|---|
| 262 | ! H2O ice |
|---|
| 263 | !-------- |
|---|
| 264 | call get_var3("h2o_ice_s",var3_read_1) |
|---|
| 265 | where (var3_read_1 < 0.) var3_read_1 = 0. |
|---|
| 266 | write(*,*) "Data for h2o_ice_s downloaded." |
|---|
| 267 | #ifndef CPP_STD |
|---|
| 268 | call get_var3("watercap",var3_read_2) |
|---|
| 269 | write(*,*) "Data for watercap downloaded." |
|---|
| 270 | #endif |
|---|
| 271 | |
|---|
| 272 | ! Compute the minimum over the year for each point |
|---|
| 273 | var2_read = minval(var3_read_1 + var3_read_2,3) |
|---|
| 274 | #ifndef CPP_1D |
|---|
| 275 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_h2o_ice(:,1,1)) |
|---|
| 276 | #else |
|---|
| 277 | min_h2o_ice(1,1,1) = var2_read(1,1) |
|---|
| 278 | #endif |
|---|
| 279 | write(*,*) "Min of h2o_ice computed." |
|---|
| 280 | |
|---|
| 281 | ! Tsurf |
|---|
| 282 | !------ |
|---|
| 283 | call get_var3("tsurf",var3_read_1) |
|---|
| 284 | write(*,*) "Data for tsurf downloaded." |
|---|
| 285 | |
|---|
| 286 | ! Compute average over the year for each point |
|---|
| 287 | var2_read = sum(var3_read_1,3)/timelen |
|---|
| 288 | #ifndef CPP_1D |
|---|
| 289 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,tsurf_avg_yr1(:,1)) |
|---|
| 290 | #else |
|---|
| 291 | tsurf_avg_yr1(1,1) = var2_read(1,1) |
|---|
| 292 | #endif |
|---|
| 293 | write(*,*) "Average of tsurf computed." |
|---|
| 294 | |
|---|
| 295 | if (soil_pem) then |
|---|
| 296 | ! Tsoil |
|---|
| 297 | !------ |
|---|
| 298 | call get_var4("soiltemp",var4_read) |
|---|
| 299 | write(*,*) "Data for soiltemp downloaded." |
|---|
| 300 | |
|---|
| 301 | ! Compute average over the year for each point |
|---|
| 302 | var3_read_3 = sum(var4_read,4)/timelen |
|---|
| 303 | #ifndef CPP_1D |
|---|
| 304 | do l = 1,nsoilmx |
|---|
| 305 | call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,var4_read(:,:,l,:),tsoil_timeseries(:,l,1,:)) |
|---|
| 306 | enddo |
|---|
| 307 | call gr_dyn_fi(nsoilmx,iim_input + 1,jjm_input + 1,ngrid,var3_read_3,tsoil_avg(:,:,1)) |
|---|
| 308 | #else |
|---|
| 309 | tsoil_timeseries(1,:,1,:) = var4_read(1,1,:,:) |
|---|
| 310 | tsoil_avg(1,:,1) = var3_read_3(1,:,1) |
|---|
| 311 | #endif |
|---|
| 312 | write(*,*) "Average of tsoil computed." |
|---|
| 313 | |
|---|
| 314 | ! Soil water density |
|---|
| 315 | !------------------- |
|---|
| 316 | call get_var4("waterdensity_soil",var4_read) |
|---|
| 317 | #ifndef CPP_1D |
|---|
| 318 | do l = 1,nsoilmx |
|---|
| 319 | call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,var4_read(:,:,l,:),watersoil_density_timeseries(:,l,1,:)) |
|---|
| 320 | enddo |
|---|
| 321 | #else |
|---|
| 322 | watersoil_density_timeseries(1,:,1,:) = var4_read(1,1,:,:) |
|---|
| 323 | #endif |
|---|
| 324 | write(*,*) "Data for waterdensity_soil downloaded." |
|---|
| 325 | |
|---|
| 326 | ! Surface water density |
|---|
| 327 | !---------------------- |
|---|
| 328 | call get_var3("waterdensity_surface",var3_read_1) |
|---|
| 329 | write(*,*) "Data for waterdensity_surface downloaded." |
|---|
| 330 | |
|---|
| 331 | ! Compute average over the year for each point |
|---|
| 332 | var2_read = sum(var3_read_1,3)/timelen |
|---|
| 333 | #ifndef CPP_1D |
|---|
| 334 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,watersurf_density_avg(:,1)) |
|---|
| 335 | #else |
|---|
| 336 | watersurf_density_avg(1,1) = var2_read(1,1) |
|---|
| 337 | #endif |
|---|
| 338 | write(*,*) "Average of waterdensity_surface computed." |
|---|
| 339 | endif |
|---|
| 340 | else ! We use slopes |
|---|
| 341 | do islope = 1,nslope |
|---|
| 342 | write(num,'(i2.2)') islope |
|---|
| 343 | |
|---|
| 344 | ! CO2 ice |
|---|
| 345 | !-------- |
|---|
| 346 | call get_var3("co2ice_slope"//num,var3_read_1) |
|---|
| 347 | where (var3_read_1 < 0.) var3_read_1 = 0. |
|---|
| 348 | write(*,*) "Data for co2_ice_slope"//num//" downloaded." |
|---|
| 349 | #ifndef CPP_STD |
|---|
| 350 | call get_var3("perennial_co2ice_slope"//num,var3_read_2) |
|---|
| 351 | write(*,*) "Data for perennial_co2ice_slope"//num//" downloaded." |
|---|
| 352 | #endif |
|---|
| 353 | |
|---|
| 354 | ! Compute the minimum over the year for each point |
|---|
| 355 | var2_read = minval(var3_read_1 + var3_read_2,3) |
|---|
| 356 | #ifndef CPP_1D |
|---|
| 357 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_co2_ice(:,islope,2)) |
|---|
| 358 | #else |
|---|
| 359 | min_co2_ice(1,islope,2) = var2_read(1,1) |
|---|
| 360 | #endif |
|---|
| 361 | write(*,*) "Min of co2_ice computed for slope "//num//"." |
|---|
| 362 | |
|---|
| 363 | ! H2O ice |
|---|
| 364 | !-------- |
|---|
| 365 | call get_var3("h2o_ice_s_slope"//num,var3_read_1) |
|---|
| 366 | where (var3_read_1 < 0.) var3_read_1 = 0. |
|---|
| 367 | write(*,*) "Data for h2o_ice_s_slope"//num//" downloaded." |
|---|
| 368 | #ifndef CPP_STD |
|---|
| 369 | call get_var3("watercap_slope"//num,var3_read_2) |
|---|
| 370 | write(*,*) "Data for watercap_slope"//num//" downloaded." |
|---|
| 371 | #endif |
|---|
| 372 | |
|---|
| 373 | ! Compute the minimum over the year for each point |
|---|
| 374 | var2_read = minval(var3_read_1 + var3_read_2,3) |
|---|
| 375 | #ifndef CPP_1D |
|---|
| 376 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_h2o_ice(:,islope,2)) |
|---|
| 377 | #else |
|---|
| 378 | min_h2o_ice(1,islope,2) = var2_read(1,1) |
|---|
| 379 | #endif |
|---|
| 380 | write(*,*) "Min of h2o_ice computed for slope "//num//"." |
|---|
| 381 | |
|---|
| 382 | ! Tsurf |
|---|
| 383 | !------ |
|---|
| 384 | call get_var3("tsurf_slope"//num,var3_read_1) |
|---|
| 385 | write(*,*) "Data for tsurf_slope"//num//" downloaded." |
|---|
| 386 | |
|---|
| 387 | ! Compute average over the year for each point |
|---|
| 388 | var2_read = sum(var3_read_1,3)/timelen |
|---|
| 389 | #ifndef CPP_1D |
|---|
| 390 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,tsurf_avg(:,islope)) |
|---|
| 391 | #else |
|---|
| 392 | tsurf_avg(1,islope) = var2_read(1,1) |
|---|
| 393 | #endif |
|---|
| 394 | write(*,*) "Average of tsurf computed for slope "//num//"." |
|---|
| 395 | |
|---|
| 396 | if (soil_pem) then |
|---|
| 397 | ! Tsoil |
|---|
| 398 | !------ |
|---|
| 399 | call get_var4("soiltemp_slope"//num,var4_read) |
|---|
| 400 | write(*,*) "Data for soiltemp_slope"//num//" downloaded." |
|---|
| 401 | |
|---|
| 402 | ! Compute average over the year for each point |
|---|
| 403 | var3_read_3 = sum(var4_read,4)/timelen |
|---|
| 404 | #ifndef CPP_1D |
|---|
| 405 | do l = 1,nsoilmx |
|---|
| 406 | call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,var4_read(:,:,l,:),tsoil_timeseries(:,l,islope,:)) |
|---|
| 407 | enddo |
|---|
| 408 | call gr_dyn_fi(nsoilmx,iim_input + 1,jjm_input + 1,ngrid,var3_read_3,tsoil_avg(:,:,islope)) |
|---|
| 409 | #else |
|---|
| 410 | tsoil_timeseries(1,:,islope,:) = var4_read(1,1,:,:) |
|---|
| 411 | tsoil_avg(1,:,islope) = var3_read_3(1,:,1) |
|---|
| 412 | #endif |
|---|
| 413 | write(*,*) "Average of tsoil computed for slope "//num//"." |
|---|
| 414 | |
|---|
| 415 | ! Soil water density |
|---|
| 416 | !------------------- |
|---|
| 417 | call get_var4("waterdensity_soil_slope"//num,var4_read) |
|---|
| 418 | #ifndef CPP_1D |
|---|
| 419 | do l = 1,nsoilmx |
|---|
| 420 | call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,var4_read(:,:,l,:),watersoil_density_timeseries(:,l,islope,:)) |
|---|
| 421 | enddo |
|---|
| 422 | #else |
|---|
| 423 | watersoil_density_timeseries(1,:,islope,:) = var4_read(1,1,:,:) |
|---|
| 424 | #endif |
|---|
| 425 | write(*,*) "Data for waterdensity_soil_slope"//num//" downloaded." |
|---|
| 426 | |
|---|
| 427 | ! Surface water density |
|---|
| 428 | !---------------------- |
|---|
| 429 | call get_var3("waterdensity_surface"//num,var3_read_1) |
|---|
| 430 | write(*,*) "Data for waterdensity_surface"//num//" downloaded." |
|---|
| 431 | |
|---|
| 432 | ! Compute average over the year for each point |
|---|
| 433 | var2_read = sum(var3_read_1,3)/timelen |
|---|
| 434 | #ifndef CPP_1D |
|---|
| 435 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,watersurf_density_avg(:,islope)) |
|---|
| 436 | #else |
|---|
| 437 | watersurf_density_avg(1,islope) = var2_read(1,1) |
|---|
| 438 | #endif |
|---|
| 439 | write(*,*) "Average of waterdensity_surface computed for slope "//num//"." |
|---|
| 440 | endif |
|---|
| 441 | enddo |
|---|
| 442 | endif |
|---|
| 443 | deallocate(var2_read,var3_read_1,var3_read_2,var3_read_3,var4_read) |
|---|
| 444 | |
|---|
| 445 | ! Close the NetCDF file |
|---|
| 446 | call error_msg(nf90_close(fID),"close",filename2) |
|---|
| 447 | write(*,*) '"> '//filename2//'" processed!' |
|---|
| 448 | |
|---|
| 449 | END SUBROUTINE read_data_PCM |
|---|
| 450 | |
|---|
| 451 | !======================================================================= |
|---|
| 452 | |
|---|
| 453 | SUBROUTINE check_dim(n1,n2,str1,str2) |
|---|
| 454 | |
|---|
| 455 | implicit none |
|---|
| 456 | |
|---|
| 457 | integer, intent(in) :: n1, n2 |
|---|
| 458 | character(len = *), intent(in) :: str1, str2 |
|---|
| 459 | |
|---|
| 460 | character(256) :: s1, s2 |
|---|
| 461 | |
|---|
| 462 | if (n1 /= n2) then |
|---|
| 463 | s1 = 'value of '//trim(str1)//' =' |
|---|
| 464 | s2 = ' read in starting file differs from parametrized '//trim(str2)//' =' |
|---|
| 465 | write(msg,'(10x,a,i4,2x,a,i4)')trim(s1),n1,trim(s2),n2 |
|---|
| 466 | call abort_gcm(trim(modname),trim(msg),1) |
|---|
| 467 | endif |
|---|
| 468 | |
|---|
| 469 | END SUBROUTINE check_dim |
|---|
| 470 | |
|---|
| 471 | !======================================================================= |
|---|
| 472 | |
|---|
| 473 | SUBROUTINE get_var1(var,v) |
|---|
| 474 | |
|---|
| 475 | implicit none |
|---|
| 476 | |
|---|
| 477 | character(len = *), intent(in) :: var |
|---|
| 478 | real, dimension(:), intent(out) :: v |
|---|
| 479 | |
|---|
| 480 | call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var) |
|---|
| 481 | call error_msg(NF90_GET_VAR(fID,vID,v),"get",var) |
|---|
| 482 | |
|---|
| 483 | END SUBROUTINE get_var1 |
|---|
| 484 | |
|---|
| 485 | !======================================================================= |
|---|
| 486 | |
|---|
| 487 | SUBROUTINE get_var3(var,v) ! on U grid |
|---|
| 488 | |
|---|
| 489 | implicit none |
|---|
| 490 | |
|---|
| 491 | character(len = *), intent(in) :: var |
|---|
| 492 | real, dimension(:,:,:), intent(out) :: v |
|---|
| 493 | |
|---|
| 494 | call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var) |
|---|
| 495 | call error_msg(NF90_GET_VAR(fID,vID,v),"get",var) |
|---|
| 496 | |
|---|
| 497 | END SUBROUTINE get_var3 |
|---|
| 498 | |
|---|
| 499 | !======================================================================= |
|---|
| 500 | |
|---|
| 501 | SUBROUTINE get_var4(var,v) |
|---|
| 502 | |
|---|
| 503 | implicit none |
|---|
| 504 | |
|---|
| 505 | character(len = *), intent(in) :: var |
|---|
| 506 | real, dimension(:,:,:,:), intent(out) :: v |
|---|
| 507 | |
|---|
| 508 | call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var) |
|---|
| 509 | call error_msg(NF90_GET_VAR(fID,vID,v),"get",var) |
|---|
| 510 | |
|---|
| 511 | END SUBROUTINE get_var4 |
|---|
| 512 | |
|---|
| 513 | !======================================================================= |
|---|
| 514 | |
|---|
| 515 | SUBROUTINE error_msg(ierr,typ,nam) |
|---|
| 516 | |
|---|
| 517 | implicit none |
|---|
| 518 | |
|---|
| 519 | integer, intent(in) :: ierr !--- NetCDF ERROR CODE |
|---|
| 520 | character(len = *), intent(in) :: typ !--- TYPE OF OPERATION |
|---|
| 521 | character(len = *), intent(in) :: nam !--- FIELD/FILE NAME |
|---|
| 522 | |
|---|
| 523 | if (ierr == nf90_noerr) return |
|---|
| 524 | select case(typ) |
|---|
| 525 | case('inq'); msg="Field <"//trim(nam)//"> is missing" |
|---|
| 526 | case('get'); msg="Reading failed for <"//trim(nam)//">" |
|---|
| 527 | case('put'); msg="Writing failed for <"//trim(nam)//">" |
|---|
| 528 | case('open'); msg="File opening failed for <"//trim(nam)//">" |
|---|
| 529 | case('close'); msg="File closing failed for <"//trim(nam)//">" |
|---|
| 530 | case default |
|---|
| 531 | write(*,*) 'There is no message for this error.' |
|---|
| 532 | error stop |
|---|
| 533 | end select |
|---|
| 534 | call abort_gcm(trim(modname),trim(msg),ierr) |
|---|
| 535 | |
|---|
| 536 | END SUBROUTINE error_msg |
|---|
| 537 | |
|---|
| 538 | END MODULE read_data_PCM_mod |
|---|
| 539 | |
|---|