m_get_mod_grid.F90 7.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169
  1. module m_get_mod_grid
  2. ! The reading of depths is not needed for assimilation of sea-ice data
  3. ! but left here for potential future use. And to retain the calling format.
  4. use netcdf
  5. #if defined (QMPI)
  6. use qmpi
  7. #else
  8. use qmpi_fake
  9. #endif
  10. private handle_err
  11. contains
  12. subroutine get_mod_grid(modlon,modlat,depths,mindx,meandx,nx,ny)
  13. implicit none
  14. ! In/out
  15. integer, intent(in) :: nx,ny ! 182, resp. 149
  16. real, dimension(nx,ny), intent(out) :: modlon, modlat, depths
  17. real, intent(out) :: mindx, meandx
  18. ! NetCDF vars
  19. integer ncid, varID, error, ncid_mask, ncid_zgr
  20. character(len=80), parameter :: maskfile = 'mask.nc' !hc!
  21. character(len=80), parameter :: meshfile_hgr = 'mesh_hgr.nc' !hc!
  22. character(len=80), parameter :: meshfile_zgr = 'mesh_zgr.nc' !hc!
  23. logical ex
  24. ! Variables for mindx ! uncomment whatever is needed
  25. real, allocatable, dimension(:,:) :: e1t, e2t!, e1u, e2u, e1v, e2v, e1f, e2f
  26. real, allocatable, dimension(:,:) :: tmask!, umask, vmask, fmask
  27. ! check the netCDF file exists
  28. inquire(file=meshfile_hgr, exist=ex)
  29. if (.not.ex) then
  30. if (master) print *, '(get_mod_grid): file does not exist: '//trim(meshfile_hgr)
  31. call stop_mpi()
  32. end if
  33. ! open the netCDF file
  34. error = nf90_open(trim(maskfile),nf90_NoWrite,ncid_mask)
  35. if (error.ne.nf90_noerr) call handle_err(error, "opening")
  36. error = nf90_open(trim(meshfile_hgr),nf90_NoWrite,ncid)
  37. if (error.ne.nf90_noerr) call handle_err(error, "opening")
  38. error = nf90_open(trim(meshfile_zgr),nf90_NoWrite,ncid_zgr)
  39. if (error.ne.nf90_noerr) call handle_err(error, "opening")
  40. ! Longitude
  41. ! Find VarID
  42. error = nf90_inq_varid(ncid, 'nav_lon', varID)
  43. if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID modlon")
  44. ! Get values of variable
  45. error = nf90_get_var(ncid, varID, modlon)
  46. if (error.ne.nf90_noerr) call handle_err(error, "getting variable modlon")
  47. ! Latitude
  48. ! Find VarID
  49. error = nf90_inq_varid(ncid, 'nav_lat', varID)
  50. if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID modlat")
  51. ! Get values of variable
  52. error = nf90_get_var(ncid, varID, modlat)
  53. if (error.ne.nf90_noerr) call handle_err(error, "getting variable modlat")
  54. ! Depths:
  55. ! Find VarID
  56. error = nf90_inq_varid(ncid_zgr, 'gdept_0', varID)
  57. if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID depths")
  58. ! Get values of variable
  59. error = nf90_get_var(ncid_zgr, varID, depths)
  60. if (error.ne.nf90_noerr) call handle_err(error, "getting variable depths")
  61. ! mindx: Smallest horizontal grid spacing. Requires some 'math'.
  62. ! Load grid spacing and corresponding masks
  63. allocate( e1t(nx,ny), e2t(nx,ny) )
  64. !allocate( e1u(nx,ny), e2u(nx,ny) ) ! In case those variables are ...
  65. !allocate( e1v(nx,ny), e2v(nx,ny) ) ! ... needed, feel free ...
  66. !allocate( e1f(nx,ny), e2f(nx,ny) ) ! ... to uncomment.
  67. allocate( tmask(nx,ny) ) ! umask(nx,ny), vmask(nx,ny), fmask(nx,ny)
  68. ! Get e1t, e1u, e1v, e1f, e2t, e2u, e2v, and e2f
  69. error = nf90_inq_varid(ncid, 'e1t', varID) ; if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID")
  70. error = nf90_get_var(ncid, varID, e1t) ; if (error.ne.nf90_noerr) call handle_err(error, "getting variable")
  71. error = nf90_inq_varid(ncid, 'e2t', varID) ; if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID")
  72. error = nf90_get_var(ncid, varID, e2t) ; if (error.ne.nf90_noerr) call handle_err(error, "getting variable")
  73. !error = nf90_inq_varid(ncid, 'e1u', varID) ; if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID")
  74. !error = nf90_get_var(ncid, varID, e1u) ; if (error.ne.nf90_noerr) call handle_err(error, "getting variable")
  75. !error = nf90_inq_varid(ncid, 'e2u', varID) ; if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID")
  76. !error = nf90_get_var(ncid, varID, e2u) ; if (error.ne.nf90_noerr) call handle_err(error, "getting variable")
  77. !error = nf90_inq_varid(ncid, 'e1v', varID) ; if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID")
  78. !error = nf90_get_var(ncid, varID, e1v) ; if (error.ne.nf90_noerr) call handle_err(error, "getting variable")
  79. !error = nf90_inq_varid(ncid, 'e2v', varID) ; if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID")
  80. !error = nf90_get_var(ncid, varID, e2v) ; if (error.ne.nf90_noerr) call handle_err(error, "getting variable")
  81. !error = nf90_inq_varid(ncid, 'e1f', varID) ; if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID")
  82. !error = nf90_get_var(ncid, varID, e1f) ; if (error.ne.nf90_noerr) call handle_err(error, "getting variable")
  83. !error = nf90_inq_varid(ncid, 'e2f', varID) ; if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID")
  84. !error = nf90_get_var(ncid, varID, e2f) ; if (error.ne.nf90_noerr) call handle_err(error, "getting variable")
  85. ! Get tmask, umask, vmask, fmask !!! only first level of 3d-nc-var fits into local var. It's all we need.
  86. error = nf90_inq_varid(ncid_mask, 'tmask', varID) ; if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID")
  87. error = nf90_get_var(ncid_mask, varID, tmask) ; if (error.ne.nf90_noerr) call handle_err(error, "getting variable")
  88. !error = nf90_inq_varid(ncid, 'umask', varID) ; if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID")
  89. !error = nf90_get_var(ncid, varID, umask) ; if (error.ne.nf90_noerr) call handle_err(error, "getting variable")
  90. !error = nf90_inq_varid(ncid, 'vmask', varID) ; if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID")
  91. !error = nf90_get_var(ncid, varID, vmask) ; if (error.ne.nf90_noerr) call handle_err(error, "getting variable")
  92. !error = nf90_inq_varid(ncid, 'fmask', varID) ; if (error.ne.nf90_noerr) call handle_err(error, "inquiring varID")
  93. !error = nf90_get_var(ncid, varID, fmask) ; if (error.ne.nf90_noerr) call handle_err(error, "getting variable")
  94. ! Smart use of min/maxval
  95. ! Find absolute minimum
  96. mindx = min(minval(e1t, tmask>0.5), minval(e2t, tmask>0.5))!, &
  97. ! minval(e1u, umask>0.5), minval(e2u, umask>0.5), &
  98. ! minval(e1v, vmask>0.5), minval(e2v, vmask>0.5), &
  99. ! minval(e1f, fmask>0.5), minval(e2f, fmask>0.5) )
  100. if (master) then
  101. print *,'(get_mod_grid) MINIMUM grid size from mesh_mask : ', mindx
  102. end if
  103. ! Find mean horizontal distance
  104. meandx = (sum(e1t,mask=tmask>0.5) + sum(e2t,mask=tmask>0.5) ) &
  105. / count(tmask>0.5)
  106. if (master) then
  107. print *,'(get_mod_grid) MEAN grid size from mesh_mask: ', meandx
  108. end if
  109. ! Safety check .. inherited from KAL
  110. if (mindx<2000.) then
  111. if (master) print *,'(get_mod_grid) min grid size lower than safety threshold - fix if you want'
  112. call stop_mpi()
  113. end if
  114. ! Safety check .. This one is not that critical so the value is set high
  115. if (mindx>500000.) then
  116. if (master) print *,'(get_mod_grid) min grid size higher than safety threshold - fix if you want'
  117. call stop_mpi()
  118. end if
  119. ! Close file
  120. error = nf90_close(ncid) ! close netCDF dataset
  121. if (error.ne.nf90_noerr) call handle_err(error, "closing")
  122. error = nf90_close(ncid_mask) ! close netCDF dataset
  123. if (error.ne.nf90_noerr) call handle_err(error, "closing")
  124. error = nf90_close(ncid_zgr) ! close netCDF dataset
  125. if (error.ne.nf90_noerr) call handle_err(error, "closing")
  126. end subroutine get_mod_grid
  127. subroutine handle_err(status, infomsg)
  128. integer, intent ( in) :: status
  129. character(len = *), intent ( in), optional :: infomsg
  130. if(status /= nf90_noerr) then
  131. if (master) then
  132. if (present(infomsg)) then
  133. print *, 'Error while '//infomsg//' - '//trim(nf90_strerror(status))
  134. else
  135. print *, trim(nf90_strerror(status))
  136. endif ! opt arg
  137. print *,'(get_mod_grid)'
  138. endif ! only master outputs
  139. call stop_mpi()
  140. end if ! check error status
  141. end subroutine handle_err
  142. end module m_get_mod_grid