m_get_mod_xyz.F90 2.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899
  1. module m_get_mod_xyz
  2. ! Gets model dimensions in file './mask.nc'
  3. ! (unless another netcdf file is submitted) and returns them.
  4. ! Added by F. Massonnet to the NERSC-ENKF routines, May 2013.
  5. ! (Presumably) Coded by C. König Beatty, in 2009
  6. ! Goal is to quickly retrieve model dimensions without using parseblk
  7. use netcdf
  8. #if defined (QMPI)
  9. use qmpi
  10. #else
  11. use qmpi_fake
  12. #endif
  13. private handle_err
  14. contains
  15. subroutine get_mod_xyz(x, y, z, moddimfilein)
  16. implicit none
  17. ! In/out
  18. integer, intent(out) :: x, y, z
  19. character(len=*), intent(in), optional :: moddimfilein
  20. ! NetCDF vars
  21. integer :: ncid, dimID, error
  22. character(len=120) :: moddimfile
  23. logical ex
  24. if (present(moddimfilein)) then
  25. moddimfile=moddimfilein
  26. else
  27. moddimfile='./mask.nc'
  28. end if
  29. ! check the netCDF file exists
  30. inquire(file=moddimfile, exist=ex)
  31. if (.not.ex) then
  32. if (master) then
  33. print *, '(get_mod_xyz): file does not exist: '//trim(moddimfile)
  34. end if
  35. call stop_mpi()
  36. end if
  37. ! open the netCDF file
  38. error = nf90_open(trim(moddimfile),nf90_NoWrite,ncid)
  39. if (error.ne.nf90_noerr) call handle_err(error, "opening")
  40. ! Find DimID of x
  41. error = nf90_inq_dimid(ncid, 'x', dimID)
  42. if (error.ne.nf90_noerr) call handle_err(error, "inquiring dimID x")
  43. ! Get size of dimension
  44. error = nf90_inquire_dimension(ncid, dimID, len = x)
  45. if (error.ne.nf90_noerr) call handle_err(error, "getting dimension x")
  46. ! Find DimID of y
  47. error = nf90_inq_dimid(ncid, 'y', dimID)
  48. if (error.ne.nf90_noerr) call handle_err(error, "inquiring dimID y")
  49. ! Get size of dimension
  50. error = nf90_inquire_dimension(ncid, dimID, len = y)
  51. if (error.ne.nf90_noerr) call handle_err(error, "getting dimension y")
  52. ! Find DimID of z
  53. error = nf90_inq_dimid(ncid, 'z', dimID)
  54. if (error.ne.nf90_noerr) call handle_err(error, "inquiring dimID z")
  55. ! Get size of dimension
  56. error = nf90_inquire_dimension(ncid, dimID, len = z)
  57. if (error.ne.nf90_noerr) call handle_err(error, "getting dimension z")
  58. ! Close file
  59. error = nf90_close(ncid) ! close netCDF dataset
  60. if (error.ne.nf90_noerr) call handle_err(error, "closing")
  61. contains
  62. subroutine handle_err(status, infomsg)
  63. integer, intent ( in) :: status
  64. character(len = *), intent ( in), optional :: infomsg
  65. if(status /= nf90_noerr) then
  66. if (present(infomsg)) then
  67. print *, 'Error while '//infomsg//' - '//trim(nf90_strerror(status))
  68. else
  69. print *, trim(nf90_strerror(status))
  70. endif
  71. stop " Stopped"
  72. end if
  73. end subroutine handle_err
  74. end subroutine get_mod_xyz
  75. end module m_get_mod_xyz