read_input_file.save.fr 8.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272
  1. module read_input_file
  2. !-----------------------------------------------------------------------
  3. !
  4. ! contains routines for reading a file to interpolate
  5. !
  6. !-----------------------------------------------------------------------
  7. use kinds_mod ! defines common data types
  8. use constants ! defines useful constants
  9. use grids ! includes all grid information
  10. use netcdf_mod ! module with netcdf vars and utilities
  11. use remap_vars ! module for all required remapping variables
  12. implicit none
  13. real (kind=dbl_kind), dimension(:,:), allocatable ::
  14. & grid1_array
  15. integer (kind=int_kind) ::
  16. & ntime, nlat, nlon, natts, nvaratts, nglobatts
  17. integer (kind=int_kind) ::
  18. & nc_infile_id, nc_time_id, nc_invar_id, dimtime
  19. real (kind=dbl_kind), dimension(:), allocatable ::
  20. & time
  21. logical :: invertlat
  22. contains
  23. !***********************************************************************
  24. subroutine read_input(infile, varname)
  25. !-----------------------------------------------------------------------
  26. !
  27. ! input variables
  28. !
  29. !-----------------------------------------------------------------------
  30. character(char_len), intent(in) ::
  31. & infile, varname ! input filename an varname
  32. !-----------------------------------------------------------------------
  33. character(char_len) ::
  34. & timename, name2
  35. real (kind=dbl_kind), dimension(:,:,:), allocatable ::
  36. & grid1_array_tmp
  37. integer, dimension(:,:,:), allocatable ::
  38. & grid1_array_tmp1
  39. real, dimension(:,:,:), allocatable ::
  40. & grid1_array_tmp2
  41. double precision, dimension(:,:,:), allocatable ::
  42. & grid1_array_tmp3
  43. integer, dimension(:), allocatable ::
  44. & time1
  45. real, dimension(:), allocatable ::
  46. & time2
  47. double precision, dimension(:), allocatable ::
  48. & time3
  49. integer (kind=int_kind) ::
  50. & ncstat, nc_var_type, ndims
  51. integer (kind=int_kind), dimension(:), allocatable :: ! netCDF ids
  52. & nc_dims_ids
  53. !-----------------------------------------------------------------------
  54. !
  55. ! read dimension information
  56. !
  57. !-----------------------------------------------------------------------
  58. ncstat = nf_open(infile, NF_NOWRITE, nc_infile_id)
  59. call netcdf_error_handler(ncstat)
  60. ncstat = nf_inq_varid(nc_infile_id, varname, nc_invar_id)
  61. call netcdf_error_handler(ncstat)
  62. ncstat = nf_inq_varndims(nc_infile_id, nc_invar_id, ndims)
  63. call netcdf_error_handler(ncstat)
  64. allocate(nc_dims_ids(ndims))
  65. ncstat = nf_inq_vardimid(nc_infile_id, nc_invar_id, nc_dims_ids)
  66. call netcdf_error_handler(ncstat)
  67. if ( ndims < 1 .or. ndims > 3) then
  68. stop "Input files should have (lon) or (lon,lat) or (lon,time)
  69. & or (lon,lat,time) dimensions"
  70. endif
  71. ncstat = nf_inq_dimlen(nc_infile_id, nc_dims_ids(1), nlon)
  72. call netcdf_error_handler(ncstat)
  73. if (nlon .ne. grid1_dims(1)) then
  74. print*," Input file : ",nlon," longitudes"
  75. print*," Weight file : ",grid1_dims(1) ," longitudes"
  76. stop " Inconsistent number of longitudes. We stop "
  77. endif
  78. name2='none'
  79. if ( ndims >= 2 ) then
  80. ncstat = nf_inq_dimname(nc_infile_id, nc_dims_ids(2), name2)
  81. endif
  82. if ( name2 /= 'time' .and. name2 /= 'time_counter' .and.
  83. & name2 /= 't' ) then
  84. ncstat = nf_inq_dimlen(nc_infile_id, nc_dims_ids(2), nlat)
  85. call netcdf_error_handler(ncstat)
  86. if (nlat .ne. grid1_dims(2)) then
  87. print*," Input file : ",nlat," latitudes"
  88. print*," Weight file : ",grid1_dims(2)," latitudes"
  89. stop " Inconsistent number of latitudes. We stop "
  90. endif
  91. dimtime=3
  92. else
  93. nlat=1
  94. if (ndims >=2) then
  95. dimtime=2
  96. endif
  97. endif
  98. !-----------------------------------------------------------------------
  99. !
  100. ! read time dimension
  101. !
  102. !-----------------------------------------------------------------------
  103. nc_time_id=-1
  104. if ( ndims == 3 .or. dimtime == 2 ) then
  105. ncstat = nf_inq_dim(nc_infile_id, nc_dims_ids(dimtime),
  106. & timename, ntime)
  107. call netcdf_error_handler(ncstat)
  108. allocate(time(ntime))
  109. allocate(time1(ntime))
  110. allocate(time2(ntime))
  111. allocate(time3(ntime))
  112. ncstat = nf_inq_varid(nc_infile_id, timename, nc_time_id)
  113. call netcdf_error_handler(ncstat)
  114. ncstat = nf_inq_vartype(nc_infile_id, nc_time_id, nc_var_type)
  115. call netcdf_error_handler(ncstat)
  116. if ( nc_var_type == NF_INT ) then
  117. ncstat = nf_get_var_int(nc_infile_id, nc_time_id, time1)
  118. if ( ncstat == NF_NOERR ) then
  119. time = real(time1,kind=dbl_kind)
  120. else
  121. stop "Problem reading input data"
  122. endif
  123. elseif (nc_var_type == NF_REAL) then
  124. ncstat = nf_get_var_real(nc_infile_id, nc_time_id, time2)
  125. if ( ncstat == NF_NOERR ) then
  126. time = real(time2,kind=dbl_kind)
  127. else
  128. stop "Problem reading input data"
  129. endif
  130. elseif ( nc_var_type == NF_DOUBLE) then
  131. ncstat = nf_get_var_double(nc_infile_id,nc_time_id, time3)
  132. if ( ncstat == NF_NOERR ) then
  133. time = real(time3,kind=dbl_kind)
  134. else
  135. stop"Problem reading input data"
  136. endif
  137. else
  138. stop"Problem with input data type"
  139. endif
  140. deallocate(time1, time2, time3)
  141. ncstat = nf_inq_varnatts(nc_infile_id, nc_time_id, natts)
  142. call netcdf_error_handler(ncstat)
  143. else
  144. ntime=1
  145. allocate(time(1))
  146. time=1.
  147. natts=0
  148. endif
  149. !-----------------------------------------------------------------------
  150. !
  151. ! allocate arrays
  152. !
  153. !-----------------------------------------------------------------------
  154. allocate( grid1_array_tmp (nlon,nlat,ntime))
  155. allocate( grid1_array_tmp1 (nlon,nlat,ntime))
  156. allocate( grid1_array_tmp2 (nlon,nlat,ntime))
  157. allocate( grid1_array_tmp3 (nlon,nlat,ntime))
  158. !-----------------------------------------------------------------------
  159. !
  160. ! read variable
  161. !
  162. !-----------------------------------------------------------------------
  163. ncstat = nf_inq_vartype(nc_infile_id, nc_invar_id, nc_var_type)
  164. call netcdf_error_handler(ncstat)
  165. if ( nc_var_type == NF_INT ) then
  166. ncstat = nf_get_var_int(nc_infile_id, nc_invar_id,
  167. & grid1_array_tmp1)
  168. if ( ncstat == NF_NOERR ) then
  169. grid1_array_tmp = real(grid1_array_tmp1,kind=dbl_kind)
  170. else
  171. stop "Problem reading input data"
  172. endif
  173. elseif (nc_var_type == NF_REAL) then
  174. ncstat = nf_get_var_real(nc_infile_id, nc_invar_id,
  175. & grid1_array_tmp2)
  176. if ( ncstat == NF_NOERR ) then
  177. grid1_array_tmp = real(grid1_array_tmp2,kind=dbl_kind)
  178. else
  179. stop "Problem reading input data"
  180. endif
  181. elseif ( nc_var_type == NF_DOUBLE) then
  182. ncstat = nf_get_var_double(nc_infile_id,nc_invar_id,
  183. & grid1_array_tmp3)
  184. if ( ncstat == NF_NOERR ) then
  185. grid1_array_tmp = real(grid1_array_tmp3,kind=dbl_kind)
  186. else
  187. stop"Problem reading input data"
  188. endif
  189. else
  190. stop"Problem with input data type"
  191. endif
  192. ncstat = nf_inq_varnatts(nc_infile_id, nc_invar_id, nvaratts)
  193. call netcdf_error_handler(ncstat)
  194. deallocate(grid1_array_tmp1, grid1_array_tmp2, grid1_array_tmp3)
  195. !-----------------------------------------------------------------------
  196. !
  197. ! reshape input file
  198. !
  199. !-----------------------------------------------------------------------
  200. allocate( grid1_array (nlon*nlat,ntime))
  201. if (invertlat) then
  202. grid1_array = reshape ( grid1_array_tmp (:,nlat:1:-1,:),
  203. & (/ nlat*nlon, ntime /) )
  204. else
  205. grid1_array = reshape (grid1_array_tmp, (/ nlat*nlon, ntime /) )
  206. endif
  207. deallocate(grid1_array_tmp)
  208. ncstat = nf_inq_varnatts(nc_infile_id, NF_GLOBAL, nglobatts)
  209. call netcdf_error_handler(ncstat)
  210. !-----------------------------------------------------------------------
  211. end subroutine read_input
  212. !***********************************************************************
  213. end module read_input_file
  214. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!