m_io_mod_fld.F90 7.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195
  1. module m_io_mod_fld
  2. ! Get or put one of the fields of a restart file, specified by
  3. ! ensemble number, field name and type, and vertical level. The
  4. ! time level is currently not used (restart files have only one)
  5. ! but who knows. Grid dimension is also needed, as well as if you
  6. ! want to 'get' or 'put'.
  7. !
  8. ! This replaces the two routines 'm_get_mod_fld.F90' and m_put_mod_fld.F90'.
  9. ! There was so much overlap that it became easier to merge the two. I think.
  10. !
  11. ! (c) July 2009, Christof.KonigBeatty@uclouvain.be
  12. use netcdf
  13. #if defined (QMPI)
  14. use qmpi
  15. #else
  16. use qmpi_fake
  17. #endif
  18. private handle_err
  19. contains
  20. subroutine io_mod_fld(fld,k,enslist,cfld,type,vlevel,tlevel,nx,ny,gorp,rdate_obs)
  21. implicit none
  22. ! In/out
  23. real,dimension(nx,ny),intent(inout):: fld ! output fl
  24. integer, intent(in) :: k ! Index to enslist
  25. integer,dimension(:), intent(in) :: enslist! List of existing ensemble members
  26. character(len=*), intent(in) :: cfld ! name of fld
  27. integer, intent(in) :: type ! which file to use
  28. integer, intent(in) :: vlevel ! vertical level (ignored)
  29. integer, intent(in) :: tlevel ! time level (ignored)
  30. integer, intent(in) :: nx,ny ! Grid dimension
  31. character(len=3), intent(in) :: gorp ! 'get' or 'put' (sorry, couldn't come up with anything better)
  32. real(kind=8), intent(in) :: rdate_obs
  33. ! NetCDF vars
  34. integer :: iens ! Ensemble member to read
  35. character(len=99) :: fcfile, anafile, cfile
  36. integer :: ncid, varID, error
  37. logical :: exfc, exan
  38. ! Other
  39. character(len=3) :: cmem
  40. integer :: zvlevel ! for i/o ocean variable
  41. ! [AD] and ice category as a dimension in SI3 2024
  42. real(kind=8) :: rdate_mod
  43. ! Find iens withing enslist
  44. iens = enslist(k)
  45. ! Create filename dep. on type of variable/parameter requested
  46. write(cmem,'(i3.3)') 100+iens ! iens=1 gives cmem = 101
  47. select case( type )
  48. case(1) ! ice variable [AD] Case with ice category as dimension (considered as "vertical" level here)
  49. fcfile ='forecast_ice_'//cmem//'.nc'
  50. anafile='analysis_ice_'//cmem//'.nc'
  51. if (vlevel>0) then
  52. zvlevel = max(vlevel,1) ! consider 3rd dimension
  53. else
  54. zvlevel = 0 ! No third dimension
  55. endif
  56. case(2) ! ocean variable
  57. fcfile ='forecast_oce_'//cmem//'.nc'
  58. anafile='analysis_oce_'//cmem//'.nc'
  59. zvlevel = max(vlevel,1)
  60. case(3) ! ice namelist parameter
  61. if (master) print *, '(io_mod_fld): ice parameter writing not implemented yet!'
  62. call stop_mpi()
  63. case(4) ! ocean namelist parameter
  64. if (master) print *, '(io_mod_fld): ocean parameter writing not implemented yet!'
  65. call stop_mpi()
  66. case default
  67. if (master) print *, '(io_mod_fld): variable type not understood!'
  68. call stop_mpi()
  69. end select
  70. ! If the fc file exists we turn it into the analysis file (unless that's already there).
  71. inquire(file=fcfile, exist=exfc)
  72. inquire(file=anafile, exist=exan)
  73. if ((.not.exfc).and.(.not.exan)) then ! Neither file is there
  74. if (master) print *, '(io_mod_fld): Restart file '//cmem//' missing!'
  75. call stop_mpi()
  76. elseif (exfc.and.(.not.exan)) then ! fcfile here but no anafile
  77. ! call system('mv '//trim(fcfile)//' '//trim(anafile) ) ! "operational" to save space
  78. call system('cp '//trim(fcfile)//' '//trim(anafile) ) ! for debugging
  79. end if
  80. ! Decide on which file to use
  81. if (gorp=='get') cfile=fcfile
  82. if (gorp=='put') cfile=anafile
  83. ! ckb prefers only one file at the time, so take care of this special case
  84. inquire(file=fcfile, exist=exfc)
  85. if (.not.exfc) cfile=anafile
  86. !!$ !XXX:
  87. !!$ write(*,*) "XXX: "
  88. !!$ write(*,*) "XXX: iens : ", iens
  89. !!$ write(*,*) "XXX: cfld : ", cfld
  90. !!$ write(*,*) "XXX: type : ", type
  91. !!$ write(*,*) "XXX: nx, ny, zvlevel: ", nx, ny, zvlevel
  92. !!$ write(*,*) "XXX: fcfile : ", trim(fcfile)
  93. !!$ write(*,*) "XXX: anafile : ", trim(anafile)
  94. !!$ write(*,*) "XXX: shape(fldIO) : ", shape(fldIO)
  95. !!$ write(*,*) "XXX: "
  96. !!$ !:XXX
  97. ! open the netCDF file
  98. error = nf90_open(trim(cfile),nf90_Write,ncid); if (error.ne.nf90_noerr) call handle_err(error, "opening")
  99. ! Find VarID of cfld
  100. error = nf90_inq_varid(ncid, trim(cfld), varID); if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID")
  101. ! Put/Get
  102. select case( type )
  103. case(3, 4) ! 2D
  104. if (gorp=='get') then
  105. error = nf90_get_var(ncid, varID, fld); if (error.ne.nf90_noerr) call handle_err(error, "getting 2D variable")
  106. elseif (gorp=='put') then
  107. error = nf90_put_var(ncid, varID, fld); if (error.ne.nf90_noerr) call handle_err(error, "putting 2D variable")
  108. else
  109. if (master) print *, "(io_mod_fld): Either 'put' or 'get'!"
  110. call stop_mpi()
  111. endif
  112. case(1) ! 2D ice varibale [AD] with possibility of 3D (ice category as a dimension // vertical dimensioni)
  113. if (zvlevel > 0) then ! 3D ice variables
  114. if (gorp=='get') then
  115. error = nf90_get_var(ncid, varID, fld,start=(/1,1,zvlevel/), count=(/nx,ny,1/))
  116. if (error.ne.nf90_noerr) call handle_err(error, "getting 3D ice variable")
  117. elseif (gorp=='put') then
  118. error = nf90_put_var(ncid, varID, fld,start=(/1,1,zvlevel/), count=(/nx,ny,1/))
  119. if (error.ne.nf90_noerr) call handle_err(error, "putting 3D ice variable")
  120. else
  121. if (master) print *, "(io_mod_fld): Either 'put' or 'get'!"
  122. call stop_mpi()
  123. endif
  124. else ! 2D ice variables
  125. if (gorp=='get') then
  126. error = nf90_get_var(ncid, varID, fld); if (error.ne.nf90_noerr) call handle_err(error, "getting 2D variable")
  127. elseif (gorp=='put') then
  128. error = nf90_put_var(ncid, varID, fld); if (error.ne.nf90_noerr) call handle_err(error, "putting 2D variable")
  129. else
  130. if (master) print *, "(io_mod_fld): Either 'put' or 'get'!"
  131. call stop_mpi()
  132. endif
  133. endif
  134. case(2) ! 3D ocean variable
  135. if (gorp=='get') then
  136. error = nf90_get_var(ncid, varID, fld, start=(/1,1,zvlevel/), count=(/nx,ny,1/))
  137. if (error.ne.nf90_noerr) call handle_err(error, "getting ocean variable")
  138. elseif (gorp=='put') then
  139. error = nf90_put_var(ncid, varID, fld, start=(/1,1,zvlevel/), count=(/nx,ny,1/))
  140. if (error.ne.nf90_noerr) call handle_err(error, "putting ocean variable")
  141. else
  142. if (master) print *, "(io_mod_fld): Either 'put' or 'get'!"
  143. call stop_mpi()
  144. endif
  145. end select
  146. !if (master) PRINT *, " Find VarID of cfld "
  147. error = nf90_inq_varid(ncid, 'time_counter', varID); if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID")
  148. error = nf90_get_var(ncid, varID, rdate_mod); if (error.ne.nf90_noerr) call handle_err(error, "getting ocean variable")
  149. ! Close file
  150. error = nf90_close(ncid); if (error.ne.nf90_noerr) call handle_err(error, "closing")
  151. ! Check date mode and date obs
  152. IF (INT(rdate_mod) .NE. INT(rdate_obs)) THEN
  153. !PRINT *, 'date mod not egal to date obs, stop, (',INT(rdate_mod),' ',INT(rdate_obs),')'
  154. !STOP 1
  155. END IF
  156. end subroutine io_mod_fld
  157. subroutine handle_err(status, infomsg)
  158. integer, intent ( in) :: status
  159. character(len = *), intent ( in), optional :: infomsg
  160. if(status /= nf90_noerr) then
  161. if (master) then
  162. if (present(infomsg)) then
  163. print *, 'Error while '//infomsg//' - '//trim(nf90_strerror(status))
  164. else
  165. print *, trim(nf90_strerror(status))
  166. endif ! opt arg
  167. print *,'(io_mod_fld)'
  168. endif ! only master outputs
  169. call stop_mpi()
  170. end if ! check error status
  171. end subroutine handle_err
  172. end module m_io_mod_fld