Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

New diagnostic grids for comparing with WOA #695

Merged
merged 2 commits into from
Aug 2, 2024
Merged
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
100 changes: 87 additions & 13 deletions src/ALE/MOM_regridding.F90
Original file line number Diff line number Diff line change
Expand Up @@ -226,12 +226,38 @@ subroutine initialize_regridding(CS, GV, US, max_depth, param_file, mdl, coord_m
real, dimension(:), allocatable :: dz_max ! Thicknesses used to find maximum interface depths
! [H ~> m or kg m-2] or other units
real, dimension(:), allocatable :: rho_target ! Target density used in HYBRID mode [kg m-3]
!> Thicknesses [m] that give level centers corresponding to table 2 of WOA09
real, dimension(40) :: woa09_dz = (/ 5., 10., 10., 15., 22.5, 25., 25., 25., &
37.5, 50., 50., 75., 100., 100., 100., 100., &
100., 100., 100., 100., 100., 100., 100., 175., &
250., 375., 500., 500., 500., 500., 500., 500., &
500., 500., 500., 500., 500., 500., 500., 500. /)
! Thicknesses [m] that give level centers approximately corresponding to table 2 of WOA09
! These are approximate because the WOA09 depths are not smoothly spaced. Levels
! 1, 4, 5, 9, 12, 24, and 36 are 2.5, 2.5, 1.25 12.5, 37.5 and 62.5 m deeper than WOA09
! but all others are identical.
real, dimension(40) :: woa09_dz_approx = (/ 5., 10., 10., 15., 22.5, 25., 25., 25., &
37.5, 50., 50., 75., 100., 100., 100., 100., &
100., 100., 100., 100., 100., 100., 100., 175., &
250., 375., 500., 500., 500., 500., 500., 500., &
500., 500., 500., 500., 500., 500., 500., 500. /)
! These are the actual spacings [m] between WOA09 depths which, if used for layer thickness, places
! the interfaces at the WOA09 depths.
real, dimension(39) :: woa09_dzi = (/ 10., 10., 10., 20., 25., 25., 25., 25., &
50., 50., 50., 100., 100., 100., 100., 100., &
100., 100., 100., 100., 100., 100., 100., 250., &
250., 500., 500., 500., 500., 500., 500., 500., &
500., 500., 500., 500., 500., 500., 500. /)
! These are the spacings [m] between WOA23 depths from table 3 of
! https://www.ncei.noaa.gov/data/oceans/woa/WOA13/DOC/woa13documentation.pdf
real, dimension(136) :: woa23_dzi = (/ 5., 5., 5., 5., 5., 5., 5., 5., 5., 5., &
5., 5., 5., 5., 5., 5., 5., 5., 5., 5., &
25., 25., 25., 25., 25., 25., 25., 25., 25., 25., &
25., 25., 25., 25., 25., 25., 50., 50., 50., 50., &
50., 50., 50., 50., 50., 50., 50., 50., 50., 50., &
50., 50., 50., 50., 50., 50., 50., 50., 50., 50., &
50., 50., 50., 50., 50., 50., 100., 100., 100., 100., &
100., 100., 100., 100., 100., 100., 100., 100., 100., 100., &
100., 100., 100., 100., 100., 100., 100., 100., 100., 100., &
100., 100., 100., 100., 100., 100., 100., 100., 100., 100., &
100., 100., 100., 100., 100., 100., 100., 100., 100., 100., &
100., 100., 100., 100., 100., 100., 100., 100., 100., 100., &
100., 100., 100., 100., 100., 100., 100., 100., 100., 100., &
100., 100., 100., 100., 100., 100. /)

call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
inputdir = slasher(inputdir)
Expand Down Expand Up @@ -325,6 +351,8 @@ subroutine initialize_regridding(CS, GV, US, max_depth, param_file, mdl, coord_m
" by a comma or space, e.g. FILE:lev.nc,dz\n"//&
" or FILE:lev.nc,interfaces=zw\n"//&
" WOA09[:N] - the WOA09 vertical grid (approximately)\n"//&
" WOA09INT[:N] - layers spanned by the WOA09 depths\n"//&
" WOA23INT[:N] - layers spanned by the WOA23 depths\n"//&
" FNC1:string - FNC1:dz_min,H_total,power,precision\n"//&
" HYBRID:string - read from a file. The string specifies\n"//&
" the filename and two variable names, separated\n"//&
Expand Down Expand Up @@ -458,29 +486,75 @@ subroutine initialize_regridding(CS, GV, US, max_depth, param_file, mdl, coord_m
call log_param(param_file, mdl, "!TARGET_DENSITIES", rho_target, &
'HYBRID target densities for interfaces', units=coordinateUnits(coord_mode))
endif
elseif (index(trim(string),'WOA09INT')==1) then
if (len_trim(string)==8) then ! string=='WOA09INT'
tmpReal = 0. ; ke = 0 ; dz_extra = 0.
do while (tmpReal<maximum_depth)
ke = ke + 1
if (ke > size(woa09_dzi)) then
dz_extra = maximum_depth - tmpReal
exit
endif
tmpReal = tmpReal + woa09_dzi(ke)
enddo
elseif (index(trim(string),'WOA09INT:')==1) then ! string starts with 'WOA09INT:'
if (len_trim(string)==9) call MOM_error(FATAL,trim(mdl)//', initialize_regridding: '// &
'Expected string of form "WOA09INT:N" but got "'//trim(string)//'".')
ke = extract_integer(string(10:len_trim(string)),'',1)
if (ke>39 .or. ke<1) call MOM_error(FATAL,trim(mdl)//', initialize_regridding: '// &
'For "WOA05INT:N" N must 0<N<40 but got "'//trim(string)//'".')
endif
allocate(dz(ke))
do k=1,min(ke, size(woa09_dzi))
dz(k) = woa09_dzi(k)
enddo
if (ke > size(woa09_dzi)) dz(ke) = dz_extra
elseif (index(trim(string),'WOA23INT')==1) then
if (len_trim(string)==8) then ! string=='WOA23INT'
tmpReal = 0. ; ke = 0 ; dz_extra = 0.
do while (tmpReal<maximum_depth)
ke = ke + 1
if (ke > size(woa23_dzi)) then
dz_extra = maximum_depth - tmpReal
exit
endif
tmpReal = tmpReal + woa23_dzi(ke)
enddo
elseif (index(trim(string),'WOA23INT:')==1) then ! string starts with 'WOA23INT:'
if (len_trim(string)==9) call MOM_error(FATAL,trim(mdl)//', initialize_regridding: '// &
'Expected string of form "WOA23INT:N" but got "'//trim(string)//'".')
ke = extract_integer(string(10:len_trim(string)),'',1)
if (ke>39 .or. ke<1) call MOM_error(FATAL,trim(mdl)//', initialize_regridding: '// &
'For "WOA05INT:N" N must 0<N<40 but got "'//trim(string)//'".')
endif
allocate(dz(ke))
do k=1,min(ke, size(woa23_dzi))
dz(k) = woa23_dzi(k)
enddo
if (ke > size(woa23_dzi)) dz(ke) = dz_extra
elseif (index(trim(string),'WOA09')==1) then
if (len_trim(string)==5) then
if (len_trim(string)==5) then ! string=='WOA09'
tmpReal = 0. ; ke = 0 ; dz_extra = 0.
do while (tmpReal<maximum_depth)
ke = ke + 1
if (ke > size(woa09_dz)) then
if (ke > size(woa09_dz_approx)) then
dz_extra = maximum_depth - tmpReal
exit
endif
tmpReal = tmpReal + woa09_dz(ke)
tmpReal = tmpReal + woa09_dz_approx(ke)
enddo
elseif (index(trim(string),'WOA09:')==1) then
elseif (index(trim(string),'WOA09:')==1) then ! string starts with 'WOA09:'
if (len_trim(string)==6) call MOM_error(FATAL,trim(mdl)//', initialize_regridding: '// &
'Expected string of form "WOA09:N" but got "'//trim(string)//'".')
ke = extract_integer(string(7:len_trim(string)),'',1)
if (ke>40 .or. ke<1) call MOM_error(FATAL,trim(mdl)//', initialize_regridding: '// &
'For "WOA05:N" N must 0<N<41 but got "'//trim(string)//'".')
endif
allocate(dz(ke))
do k=1,min(ke, size(woa09_dz))
dz(k) = woa09_dz(k)
do k=1,min(ke, size(woa09_dz_approx))
dz(k) = woa09_dz_approx(k)
enddo
if (ke > size(woa09_dz)) dz(ke) = dz_extra
if (ke > size(woa09_dz_approx)) dz(ke) = dz_extra
else
call MOM_error(FATAL,trim(mdl)//", initialize_regridding: "// &
"Unrecognized coordinate configuration"//trim(string))
Expand Down
Loading