[3096] | 1 | MODULE read_data_PCM_mod |
---|
[2794] | 2 | |
---|
[3076] | 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 |
---|
[2794] | 6 | |
---|
[3076] | 7 | implicit none |
---|
| 8 | |
---|
[3122] | 9 | character(13), parameter :: modname = 'read_data_PCM' |
---|
| 10 | character(256) :: msg ! for reading |
---|
| 11 | integer :: fID, vID ! for reading |
---|
[3076] | 12 | |
---|
[2794] | 13 | !======================================================================= |
---|
[3076] | 14 | contains |
---|
| 15 | !======================================================================= |
---|
| 16 | |
---|
[3571] | 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) |
---|
[3076] | 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 | !======================================================================= |
---|
[2794] | 26 | ! |
---|
[3096] | 27 | ! Purpose: Read initial confitions file from the PCM |
---|
[2794] | 28 | ! |
---|
[3602] | 29 | ! Authors: JBC |
---|
[2794] | 30 | !======================================================================= |
---|
| 31 | |
---|
[3076] | 32 | include "dimensions.h" |
---|
[2794] | 33 | |
---|
[3076] | 34 | !======================================================================= |
---|
[3571] | 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 |
---|
[2855] | 39 | ! Ouputs |
---|
[3571] | 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 |
---|
[3620] | 54 | integer :: i, j, t, islope ! Loop variables |
---|
[3609] | 55 | real :: A, B ! Intermediate variables to compute the mean molar mass of the layer |
---|
| 56 | character(:), allocatable :: num ! For reading sloped variables |
---|
[3602] | 57 | real, dimension(:,:), allocatable :: var2_read ! Variables for reading (2 dimensions) |
---|
| 58 | real, dimension(:,:,:), allocatable :: var3_read_1, var3_read_2, var3_read_3 ! Variables for reading (3 dimensions) |
---|
| 59 | real, dimension(:,:,:,:), allocatable :: var4_read ! Variables for reading (4 dimensions) |
---|
[2794] | 60 | !----------------------------------------------------------------------- |
---|
[3571] | 61 | |
---|
[3602] | 62 | !------------------------------- Year 1 -------------------------------- |
---|
[3122] | 63 | ! Open the NetCDF file |
---|
[3571] | 64 | write(*,*) "Opening "//filename1//"..." |
---|
| 65 | call error_msg(NF90_OPEN(filename1,NF90_NOWRITE,fID),"open",filename1) |
---|
[2794] | 66 | |
---|
[3598] | 67 | ! Download the data from the file |
---|
[3602] | 68 | 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)) |
---|
| 69 | |
---|
[3609] | 70 | if (nslope == 1) then ! No slope |
---|
| 71 | allocate(character(0) :: num) |
---|
| 72 | else ! We use slopes |
---|
| 73 | allocate(character(8) :: num) |
---|
| 74 | endif |
---|
| 75 | |
---|
| 76 | do islope = 1,nslope |
---|
| 77 | if (nslope /= 1) then |
---|
[3611] | 78 | num=' ' |
---|
[3609] | 79 | write(num,'(i2.2)') islope |
---|
| 80 | num = '_slope'//num |
---|
| 81 | endif |
---|
| 82 | |
---|
[3602] | 83 | ! CO2 ice |
---|
| 84 | !-------- |
---|
[3609] | 85 | call get_var3("co2ice"//num,var3_read_1) |
---|
[3602] | 86 | where (var3_read_1 < 0.) var3_read_1 = 0. |
---|
[3609] | 87 | write(*,*) "Data for co2_ice"//num//" downloaded." |
---|
[3602] | 88 | #ifndef CPP_STD |
---|
[3609] | 89 | call get_var3("perennial_co2ice"//num,var3_read_2) |
---|
| 90 | write(*,*) "Data for perennial_co2ice"//num//" downloaded." |
---|
[3602] | 91 | #endif |
---|
[3571] | 92 | |
---|
[3602] | 93 | ! Compute the minimum over the year for each point |
---|
| 94 | var2_read = minval(var3_read_1 + var3_read_2,3) |
---|
| 95 | #ifndef CPP_1D |
---|
[3609] | 96 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_co2_ice(:,islope,1)) |
---|
[3602] | 97 | #else |
---|
[3609] | 98 | min_co2_ice(1,islope,1) = var2_read(1,1) |
---|
[3602] | 99 | #endif |
---|
[3609] | 100 | write(*,*) "Min of co2_ice"//num//" computed." |
---|
[3602] | 101 | |
---|
| 102 | ! H2O ice |
---|
| 103 | !-------- |
---|
[3609] | 104 | call get_var3("h2o_ice_s"//num,var3_read_1) |
---|
[3602] | 105 | where (var3_read_1 < 0.) var3_read_1 = 0. |
---|
[3609] | 106 | write(*,*) "Data for h2o_ice_s"//num//" downloaded." |
---|
[3571] | 107 | #ifndef CPP_STD |
---|
[3609] | 108 | call get_var3("watercap"//num,var3_read_2) |
---|
| 109 | write(*,*) "Data for watercap"//num//" downloaded." |
---|
[3602] | 110 | #endif |
---|
[3571] | 111 | |
---|
[3602] | 112 | ! Compute the minimum over the year for each point |
---|
| 113 | var2_read = minval(var3_read_1 + var3_read_2,3) |
---|
| 114 | #ifndef CPP_1D |
---|
[3609] | 115 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_h2o_ice(:,islope,1)) |
---|
[3602] | 116 | #else |
---|
[3609] | 117 | min_h2o_ice(1,islope,1) = var2_read(1,1) |
---|
[3571] | 118 | #endif |
---|
[3609] | 119 | write(*,*) "Min of h2o_ice"//num//" computed." |
---|
[3602] | 120 | |
---|
| 121 | ! Tsurf |
---|
| 122 | !------ |
---|
[3609] | 123 | call get_var3("tsurf"//num,var3_read_1) |
---|
| 124 | write(*,*) "Data for tsurf"//num//" downloaded." |
---|
[3602] | 125 | |
---|
| 126 | ! Compute average over the year for each point |
---|
| 127 | var2_read = sum(var3_read_1,3)/timelen |
---|
| 128 | #ifndef CPP_1D |
---|
[3609] | 129 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,tsurf_avg_yr1(:,islope)) |
---|
[3602] | 130 | #else |
---|
[3609] | 131 | tsurf_avg_yr1(1,islope) = var2_read(1,1) |
---|
[3602] | 132 | #endif |
---|
[3609] | 133 | write(*,*) "Average of tsurf"//num//" computed." |
---|
| 134 | enddo |
---|
[3602] | 135 | |
---|
| 136 | ! Close the NetCDF file |
---|
| 137 | call error_msg(nf90_close(fID),"close",filename1) |
---|
[3603] | 138 | write(*,*) '> "'//filename1//'" processed!' |
---|
[3571] | 139 | |
---|
[3602] | 140 | !------------------------------- Year 2 -------------------------------- |
---|
[3571] | 141 | ! Open the NetCDF file |
---|
| 142 | write(*,*) "Opening "//filename2//"..." |
---|
| 143 | call error_msg(NF90_OPEN(filename2,NF90_NOWRITE,fID),"open",filename2) |
---|
| 144 | |
---|
[3363] | 145 | ! Download the data from the file |
---|
[3602] | 146 | allocate(var3_read_3(iim_input + 1,jjm_input + 1,nsoilmx),var4_read(iim_input + 1,jjm_input + 1,nsoilmx,timelen)) |
---|
| 147 | |
---|
| 148 | ! Surface pressure |
---|
| 149 | !----------------- |
---|
| 150 | call get_var3("ps",var3_read_1) |
---|
[3591] | 151 | write(*,*) "Data for surface pressure downloaded." |
---|
[2794] | 152 | |
---|
[3602] | 153 | ! Compute average over the year for each point |
---|
| 154 | var2_read = sum(var3_read_1,3)/timelen |
---|
| 155 | #ifndef CPP_1D |
---|
| 156 | call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,var3_read_1,ps_timeseries) |
---|
| 157 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,ps_avg) |
---|
| 158 | #else |
---|
[3609] | 159 | ps_timeseries(1,:) = var3_read_1(1,1,:) |
---|
[3602] | 160 | ps_avg(1) = var2_read(1,1) |
---|
| 161 | #endif |
---|
| 162 | write(*,*) "Average of surface pressure computed." |
---|
[2794] | 163 | |
---|
[3602] | 164 | ! CO2 vmr |
---|
| 165 | !-------- |
---|
| 166 | call get_var3("co2_layer1",var3_read_1) |
---|
| 167 | where (var3_read_1 < 0.) |
---|
| 168 | var3_read_1 = 1.e-10 |
---|
| 169 | else where (var3_read_1 > 1.) |
---|
| 170 | var3_read_1 = 1. |
---|
| 171 | end where |
---|
| 172 | #ifndef CPP_1D |
---|
| 173 | call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,var3_read_1,q_co2) |
---|
| 174 | #else |
---|
| 175 | q_co2(1,:) = var3_read_1(1,1,:) |
---|
| 176 | #endif |
---|
| 177 | A = (1./m_co2 - 1./m_noco2) |
---|
| 178 | B = 1./m_noco2 |
---|
| 179 | vmr_co2_PCM = q_co2/(A*q_co2 + B)/m_co2 |
---|
| 180 | write(*,*) "Data for CO2 vmr downloaded." |
---|
[2794] | 181 | |
---|
[3602] | 182 | ! H2O vmr |
---|
| 183 | !-------- |
---|
| 184 | call get_var3("h2o_layer1",var3_read_1) |
---|
| 185 | where (var3_read_1 < 0.) |
---|
| 186 | var3_read_1 = 1.e-10 |
---|
| 187 | else where (var3_read_1 > 1.) |
---|
| 188 | var3_read_1 = 1. |
---|
| 189 | end where |
---|
| 190 | #ifndef CPP_1D |
---|
| 191 | call gr_dyn_fi(timelen,iim_input + 1,jjm_input + 1,ngrid,var3_read_1,q_h2o) |
---|
| 192 | #else |
---|
| 193 | q_h2o(1,:) = var3_read_1(1,1,:) |
---|
| 194 | #endif |
---|
| 195 | write(*,*) "Data for H2O vmr downloaded." |
---|
| 196 | |
---|
[3609] | 197 | do islope = 1,nslope |
---|
| 198 | if (nslope /= 1) then |
---|
[3611] | 199 | num=' ' |
---|
[3609] | 200 | write(num,'(i2.2)') islope |
---|
| 201 | num = '_slope'//num |
---|
| 202 | endif |
---|
| 203 | |
---|
[3602] | 204 | ! CO2 ice |
---|
| 205 | !-------- |
---|
[3609] | 206 | call get_var3("co2ice"//num,var3_read_1) |
---|
[3602] | 207 | where (var3_read_1 < 0.) var3_read_1 = 0. |
---|
[3609] | 208 | write(*,*) "Data for co2_ice"//num//" downloaded." |
---|
[3602] | 209 | #ifndef CPP_STD |
---|
[3609] | 210 | call get_var3("perennial_co2ice"//num,var3_read_2) |
---|
| 211 | write(*,*) "Data for perennial_co2ice"//num//" downloaded." |
---|
[3602] | 212 | #endif |
---|
[2835] | 213 | |
---|
[3602] | 214 | ! Compute the minimum over the year for each point |
---|
| 215 | var2_read = minval(var3_read_1 + var3_read_2,3) |
---|
| 216 | #ifndef CPP_1D |
---|
[3609] | 217 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_co2_ice(:,islope,2)) |
---|
[3602] | 218 | #else |
---|
[3609] | 219 | min_co2_ice(1,islope,2) = var2_read(1,1) |
---|
[3602] | 220 | #endif |
---|
[3609] | 221 | write(*,*) "Min of co2_ice"//num//" computed." |
---|
[3602] | 222 | |
---|
| 223 | ! H2O ice |
---|
| 224 | !-------- |
---|
[3609] | 225 | call get_var3("h2o_ice_s"//num,var3_read_1) |
---|
[3602] | 226 | where (var3_read_1 < 0.) var3_read_1 = 0. |
---|
[3609] | 227 | write(*,*) "Data for h2o_ice_s"//num//" downloaded." |
---|
[3602] | 228 | #ifndef CPP_STD |
---|
[3609] | 229 | call get_var3("watercap"//num,var3_read_2) |
---|
| 230 | write(*,*) "Data for watercap"//num//" downloaded." |
---|
[3602] | 231 | #endif |
---|
[2835] | 232 | |
---|
[3602] | 233 | ! Compute the minimum over the year for each point |
---|
| 234 | var2_read = minval(var3_read_1 + var3_read_2,3) |
---|
| 235 | #ifndef CPP_1D |
---|
[3609] | 236 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,min_h2o_ice(:,islope,2)) |
---|
[3602] | 237 | #else |
---|
[3609] | 238 | min_h2o_ice(1,islope,2) = var2_read(1,1) |
---|
[3602] | 239 | #endif |
---|
[3609] | 240 | write(*,*) "Min of h2o_ice"//num//" computed." |
---|
[3602] | 241 | |
---|
| 242 | ! Tsurf |
---|
| 243 | !------ |
---|
[3609] | 244 | call get_var3("tsurf"//num,var3_read_1) |
---|
[3620] | 245 | write(*,*) "Data for tsurf"//num//" downloaded." |
---|
[3143] | 246 | |
---|
[3602] | 247 | ! Compute average over the year for each point |
---|
| 248 | var2_read = sum(var3_read_1,3)/timelen |
---|
| 249 | #ifndef CPP_1D |
---|
[3609] | 250 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,tsurf_avg(:,islope)) |
---|
[3602] | 251 | #else |
---|
[3609] | 252 | tsurf_avg(1,islope) = var2_read(1,1) |
---|
[3602] | 253 | #endif |
---|
[3609] | 254 | write(*,*) "Average of tsurf"//num//" computed." |
---|
[3130] | 255 | |
---|
[3122] | 256 | if (soil_pem) then |
---|
[3602] | 257 | ! Tsoil |
---|
| 258 | !------ |
---|
[3609] | 259 | call get_var4("soiltemp"//num,var4_read) |
---|
| 260 | write(*,*) "Data for soiltemp"//num//" downloaded." |
---|
[2794] | 261 | |
---|
[3602] | 262 | ! Compute average over the year for each point |
---|
| 263 | var3_read_3 = sum(var4_read,4)/timelen |
---|
| 264 | #ifndef CPP_1D |
---|
[3620] | 265 | do t = 1,timelen |
---|
| 266 | call gr_dyn_fi(nsoilmx,iim_input + 1,jjm_input + 1,ngrid,var4_read(:,:,:,t),tsoil_timeseries(:,:,islope,t)) |
---|
[3602] | 267 | enddo |
---|
[3609] | 268 | call gr_dyn_fi(nsoilmx,iim_input + 1,jjm_input + 1,ngrid,var3_read_3,tsoil_avg(:,:,islope)) |
---|
[3602] | 269 | #else |
---|
[3609] | 270 | tsoil_timeseries(1,:,islope,:) = var4_read(1,1,:,:) |
---|
| 271 | tsoil_avg(1,:,islope) = var3_read_3(1,1,:) |
---|
[3602] | 272 | #endif |
---|
[3609] | 273 | write(*,*) "Average of tsoil"//num//" computed." |
---|
[3602] | 274 | |
---|
| 275 | ! Soil water density |
---|
| 276 | !------------------- |
---|
[3609] | 277 | call get_var4("waterdensity_soil"//num,var4_read) |
---|
[3602] | 278 | #ifndef CPP_1D |
---|
[3620] | 279 | do t = 1,timelen |
---|
| 280 | call gr_dyn_fi(nsoilmx,iim_input + 1,jjm_input + 1,ngrid,var4_read(:,:,:,t),watersoil_density_timeseries(:,:,islope,t)) |
---|
[3602] | 281 | enddo |
---|
| 282 | #else |
---|
[3609] | 283 | watersoil_density_timeseries(1,:,islope,:) = var4_read(1,1,:,:) |
---|
[3602] | 284 | #endif |
---|
[3609] | 285 | write(*,*) "Data for waterdensity_soil"//num//" downloaded." |
---|
[2794] | 286 | |
---|
[3602] | 287 | ! Surface water density |
---|
| 288 | !---------------------- |
---|
[3609] | 289 | call get_var3("waterdensity_surface"//num,var3_read_1) |
---|
| 290 | write(*,*) "Data for waterdensity_surface"//num//" downloaded." |
---|
[3602] | 291 | |
---|
| 292 | ! Compute average over the year for each point |
---|
| 293 | var2_read = sum(var3_read_1,3)/timelen |
---|
| 294 | #ifndef CPP_1D |
---|
[3609] | 295 | call gr_dyn_fi(1,iim_input + 1,jjm_input + 1,ngrid,var2_read,watersurf_density_avg(:,islope)) |
---|
[3602] | 296 | #else |
---|
[3609] | 297 | watersurf_density_avg(1,islope) = var2_read(1,1) |
---|
[3122] | 298 | #endif |
---|
[3609] | 299 | write(*,*) "Average of waterdensity_surface"//num//" computed." |
---|
[3602] | 300 | endif |
---|
[3609] | 301 | enddo |
---|
| 302 | deallocate(var2_read,var3_read_1,var3_read_2,var3_read_3,var4_read,num) |
---|
[3602] | 303 | |
---|
| 304 | ! Close the NetCDF file |
---|
| 305 | call error_msg(nf90_close(fID),"close",filename2) |
---|
[3616] | 306 | write(*,*) '> "'//filename2//'" processed!' |
---|
[3602] | 307 | |
---|
[3096] | 308 | END SUBROUTINE read_data_PCM |
---|
[2897] | 309 | |
---|
[3122] | 310 | !======================================================================= |
---|
| 311 | |
---|
[3076] | 312 | SUBROUTINE check_dim(n1,n2,str1,str2) |
---|
[2980] | 313 | |
---|
[3076] | 314 | implicit none |
---|
[2980] | 315 | |
---|
[3076] | 316 | integer, intent(in) :: n1, n2 |
---|
| 317 | character(len = *), intent(in) :: str1, str2 |
---|
[2980] | 318 | |
---|
[3076] | 319 | character(256) :: s1, s2 |
---|
[2794] | 320 | |
---|
[3076] | 321 | if (n1 /= n2) then |
---|
| 322 | s1 = 'value of '//trim(str1)//' =' |
---|
| 323 | s2 = ' read in starting file differs from parametrized '//trim(str2)//' =' |
---|
| 324 | write(msg,'(10x,a,i4,2x,a,i4)')trim(s1),n1,trim(s2),n2 |
---|
| 325 | call abort_gcm(trim(modname),trim(msg),1) |
---|
| 326 | endif |
---|
| 327 | |
---|
[2794] | 328 | END SUBROUTINE check_dim |
---|
| 329 | |
---|
[3076] | 330 | !======================================================================= |
---|
[2794] | 331 | |
---|
[3130] | 332 | SUBROUTINE get_var1(var,v) |
---|
[3076] | 333 | |
---|
| 334 | implicit none |
---|
| 335 | |
---|
[3130] | 336 | character(len = *), intent(in) :: var |
---|
[3076] | 337 | real, dimension(:), intent(out) :: v |
---|
| 338 | |
---|
[3130] | 339 | call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var) |
---|
| 340 | call error_msg(NF90_GET_VAR(fID,vID,v),"get",var) |
---|
[3076] | 341 | |
---|
[2794] | 342 | END SUBROUTINE get_var1 |
---|
| 343 | |
---|
[3076] | 344 | !======================================================================= |
---|
[2794] | 345 | |
---|
[3130] | 346 | SUBROUTINE get_var3(var,v) ! on U grid |
---|
[2794] | 347 | |
---|
[3076] | 348 | implicit none |
---|
| 349 | |
---|
[3130] | 350 | character(len = *), intent(in) :: var |
---|
[3076] | 351 | real, dimension(:,:,:), intent(out) :: v |
---|
| 352 | |
---|
[3130] | 353 | call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var) |
---|
| 354 | call error_msg(NF90_GET_VAR(fID,vID,v),"get",var) |
---|
[3076] | 355 | |
---|
[2794] | 356 | END SUBROUTINE get_var3 |
---|
| 357 | |
---|
[3076] | 358 | !======================================================================= |
---|
| 359 | |
---|
[3130] | 360 | SUBROUTINE get_var4(var,v) |
---|
[3076] | 361 | |
---|
| 362 | implicit none |
---|
| 363 | |
---|
[3130] | 364 | character(len = *), intent(in) :: var |
---|
[3076] | 365 | real, dimension(:,:,:,:), intent(out) :: v |
---|
| 366 | |
---|
[3130] | 367 | call error_msg(NF90_INQ_VARID(fID,var,vID),"inq",var) |
---|
| 368 | call error_msg(NF90_GET_VAR(fID,vID,v),"get",var) |
---|
[3076] | 369 | |
---|
[2794] | 370 | END SUBROUTINE get_var4 |
---|
| 371 | |
---|
[3076] | 372 | !======================================================================= |
---|
[2794] | 373 | |
---|
[3076] | 374 | SUBROUTINE error_msg(ierr,typ,nam) |
---|
| 375 | |
---|
| 376 | implicit none |
---|
| 377 | |
---|
| 378 | integer, intent(in) :: ierr !--- NetCDF ERROR CODE |
---|
| 379 | character(len = *), intent(in) :: typ !--- TYPE OF OPERATION |
---|
| 380 | character(len = *), intent(in) :: nam !--- FIELD/FILE NAME |
---|
| 381 | |
---|
| 382 | if (ierr == nf90_noerr) return |
---|
| 383 | select case(typ) |
---|
| 384 | case('inq'); msg="Field <"//trim(nam)//"> is missing" |
---|
| 385 | case('get'); msg="Reading failed for <"//trim(nam)//">" |
---|
[3363] | 386 | case('put'); msg="Writing failed for <"//trim(nam)//">" |
---|
[3076] | 387 | case('open'); msg="File opening failed for <"//trim(nam)//">" |
---|
| 388 | case('close'); msg="File closing failed for <"//trim(nam)//">" |
---|
| 389 | case default |
---|
| 390 | write(*,*) 'There is no message for this error.' |
---|
| 391 | error stop |
---|
| 392 | end select |
---|
| 393 | call abort_gcm(trim(modname),trim(msg),ierr) |
---|
| 394 | |
---|
| 395 | END SUBROUTINE error_msg |
---|
| 396 | |
---|
[3096] | 397 | END MODULE read_data_PCM_mod |
---|
[3602] | 398 | |
---|