read_input_file.f 8.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276
  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, ji
  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. nlat=1
  80. if ( ndims >= 2 ) then
  81. ncstat = nf_inq_dimname(nc_infile_id, nc_dims_ids(2), name2)
  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. dimtime=2
  94. endif
  95. endif
  96. !-----------------------------------------------------------------------
  97. !
  98. ! read time dimension
  99. !
  100. !-----------------------------------------------------------------------
  101. nc_time_id=-1
  102. if ( ndims == 3 .or. dimtime == 2 ) then
  103. ncstat = nf_inq_dim(nc_infile_id, nc_dims_ids(dimtime),
  104. & timename, ntime)
  105. call netcdf_error_handler(ncstat)
  106. allocate(time(ntime))
  107. ncstat = nf_inq_varid(nc_infile_id, timename, nc_time_id)
  108. if (ncstat == NF_NOERR) then
  109. ncstat = nf_inq_vartype(nc_infile_id, nc_time_id, nc_var_type)
  110. call netcdf_error_handler(ncstat)
  111. allocate(time1(ntime))
  112. allocate(time2(ntime))
  113. allocate(time3(ntime))
  114. if ( nc_var_type == NF_INT ) then
  115. ncstat = nf_get_var_int(nc_infile_id, nc_time_id, time1)
  116. if ( ncstat == NF_NOERR ) then
  117. time = real(time1,kind=dbl_kind)
  118. else
  119. stop "Problem reading input data"
  120. endif
  121. elseif (nc_var_type == NF_REAL) then
  122. ncstat = nf_get_var_real(nc_infile_id, nc_time_id, time2)
  123. if ( ncstat == NF_NOERR ) then
  124. time = real(time2,kind=dbl_kind)
  125. else
  126. stop "Problem reading input data"
  127. endif
  128. elseif ( nc_var_type == NF_DOUBLE) then
  129. ncstat = nf_get_var_double(nc_infile_id,nc_time_id, time3)
  130. if ( ncstat == NF_NOERR ) then
  131. time = real(time3,kind=dbl_kind)
  132. else
  133. stop"Problem reading input data"
  134. endif
  135. else
  136. stop"Problem with input data type"
  137. endif
  138. deallocate(time1, time2, time3)
  139. ncstat = nf_inq_varnatts(nc_infile_id, nc_time_id, natts)
  140. call netcdf_error_handler(ncstat)
  141. else
  142. do ji=1,ntime
  143. time(ji)=ji
  144. enddo
  145. natts=0
  146. endif
  147. else
  148. ntime=1
  149. allocate(time(1))
  150. time=1.
  151. natts=0
  152. endif
  153. !-----------------------------------------------------------------------
  154. !
  155. ! allocate arrays
  156. !
  157. !-----------------------------------------------------------------------
  158. allocate( grid1_array_tmp (nlon,nlat,ntime))
  159. allocate( grid1_array_tmp1 (nlon,nlat,ntime))
  160. allocate( grid1_array_tmp2 (nlon,nlat,ntime))
  161. allocate( grid1_array_tmp3 (nlon,nlat,ntime))
  162. !-----------------------------------------------------------------------
  163. !
  164. ! read variable
  165. !
  166. !-----------------------------------------------------------------------
  167. ncstat = nf_inq_vartype(nc_infile_id, nc_invar_id, nc_var_type)
  168. call netcdf_error_handler(ncstat)
  169. if ( nc_var_type == NF_INT ) then
  170. ncstat = nf_get_var_int(nc_infile_id, nc_invar_id,
  171. & grid1_array_tmp1)
  172. if ( ncstat == NF_NOERR ) then
  173. grid1_array_tmp = real(grid1_array_tmp1,kind=dbl_kind)
  174. else
  175. stop "Problem reading input data"
  176. endif
  177. elseif (nc_var_type == NF_REAL) then
  178. ncstat = nf_get_var_real(nc_infile_id, nc_invar_id,
  179. & grid1_array_tmp2)
  180. if ( ncstat == NF_NOERR ) then
  181. grid1_array_tmp = real(grid1_array_tmp2,kind=dbl_kind)
  182. else
  183. stop "Problem reading input data"
  184. endif
  185. elseif ( nc_var_type == NF_DOUBLE) then
  186. ncstat = nf_get_var_double(nc_infile_id,nc_invar_id,
  187. & grid1_array_tmp3)
  188. if ( ncstat == NF_NOERR ) then
  189. grid1_array_tmp = real(grid1_array_tmp3,kind=dbl_kind)
  190. else
  191. stop"Problem reading input data"
  192. endif
  193. else
  194. stop"Problem with input data type"
  195. endif
  196. ncstat = nf_inq_varnatts(nc_infile_id, nc_invar_id, nvaratts)
  197. call netcdf_error_handler(ncstat)
  198. deallocate(grid1_array_tmp1, grid1_array_tmp2, grid1_array_tmp3)
  199. !-----------------------------------------------------------------------
  200. !
  201. ! reshape input file
  202. !
  203. !-----------------------------------------------------------------------
  204. allocate( grid1_array (nlon*nlat,ntime))
  205. if (invertlat) then
  206. grid1_array = reshape ( grid1_array_tmp (:,nlat:1:-1,:),
  207. & (/ nlat*nlon, ntime /) )
  208. else
  209. grid1_array = reshape (grid1_array_tmp, (/ nlat*nlon, ntime /) )
  210. endif
  211. deallocate(grid1_array_tmp)
  212. ncstat = nf_inq_varnatts(nc_infile_id, NF_GLOBAL, nglobatts)
  213. call netcdf_error_handler(ncstat)
  214. !-----------------------------------------------------------------------
  215. end subroutine read_input
  216. !***********************************************************************
  217. end module read_input_file
  218. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!