MOM6
coord_hycom.F90
1 !> Regrid columns for the HyCOM coordinate
2 module coord_hycom
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_error_handler, only : mom_error, fatal
8 use regrid_interp, only : interp_cs_type, build_and_interpolate_grid
9 
10 implicit none ; private
11 
12 !> Control structure containing required parameters for the HyCOM coordinate
13 type, public :: hycom_cs ; private
14 
15  !> Number of layers/levels in generated grid
16  integer :: nk
17 
18  !> Nominal near-surface resolution [Z ~> m]
19  real, allocatable, dimension(:) :: coordinateresolution
20 
21  !> Nominal density of interfaces [R ~> kg m-3]
22  real, allocatable, dimension(:) :: target_density
23 
24  !> Maximum depths of interfaces [H ~> m or kg m-2]
25  real, allocatable, dimension(:) :: max_interface_depths
26 
27  !> Maximum thicknesses of layers [H ~> m or kg m-2]
28  real, allocatable, dimension(:) :: max_layer_thickness
29 
30  !> Interpolation control structure
31  type(interp_cs_type) :: interp_cs
32 end type hycom_cs
33 
34 public init_coord_hycom, set_hycom_params, build_hycom1_column, end_coord_hycom
35 
36 contains
37 
38 !> Initialise a hycom_CS with pointers to parameters
39 subroutine init_coord_hycom(CS, nk, coordinateResolution, target_density, interp_CS)
40  type(hycom_cs), pointer :: cs !< Unassociated pointer to hold the control structure
41  integer, intent(in) :: nk !< Number of layers in generated grid
42  real, dimension(nk), intent(in) :: coordinateresolution !< Nominal near-surface resolution [Z ~> m]
43  real, dimension(nk+1),intent(in) :: target_density !< Interface target densities [R ~> kg m-3]
44  type(interp_cs_type), intent(in) :: interp_cs !< Controls for interpolation
45 
46  if (associated(cs)) call mom_error(fatal, "init_coord_hycom: CS already associated!")
47  allocate(cs)
48  allocate(cs%coordinateResolution(nk))
49  allocate(cs%target_density(nk+1))
50 
51  cs%nk = nk
52  cs%coordinateResolution(:) = coordinateresolution(:)
53  cs%target_density(:) = target_density(:)
54  cs%interp_CS = interp_cs
55 
56 end subroutine init_coord_hycom
57 
58 !> This subroutine deallocates memory in the control structure for the coord_hycom module
59 subroutine end_coord_hycom(CS)
60  type(hycom_cs), pointer :: cs !< Coordinate control structure
61 
62  ! nothing to do
63  if (.not. associated(cs)) return
64  deallocate(cs%coordinateResolution)
65  deallocate(cs%target_density)
66  if (allocated(cs%max_interface_depths)) deallocate(cs%max_interface_depths)
67  if (allocated(cs%max_layer_thickness)) deallocate(cs%max_layer_thickness)
68  deallocate(cs)
69 end subroutine end_coord_hycom
70 
71 !> This subroutine can be used to set the parameters for the coord_hycom module
72 subroutine set_hycom_params(CS, max_interface_depths, max_layer_thickness, interp_CS)
73  type(hycom_cs), pointer :: cs !< Coordinate control structure
74  real, dimension(:), optional, intent(in) :: max_interface_depths !< Maximum depths of interfaces [H ~> m or kg m-2]
75  real, dimension(:), optional, intent(in) :: max_layer_thickness !< Maximum thicknesses of layers [H ~> m or kg m-2]
76  type(interp_cs_type), optional, intent(in) :: interp_cs !< Controls for interpolation
77 
78  if (.not. associated(cs)) call mom_error(fatal, "set_hycom_params: CS not associated")
79 
80  if (present(max_interface_depths)) then
81  if (size(max_interface_depths) /= cs%nk+1) &
82  call mom_error(fatal, "set_hycom_params: max_interface_depths inconsistent size")
83  allocate(cs%max_interface_depths(cs%nk+1))
84  cs%max_interface_depths(:) = max_interface_depths(:)
85  endif
86 
87  if (present(max_layer_thickness)) then
88  if (size(max_layer_thickness) /= cs%nk) &
89  call mom_error(fatal, "set_hycom_params: max_layer_thickness inconsistent size")
90  allocate(cs%max_layer_thickness(cs%nk))
91  cs%max_layer_thickness(:) = max_layer_thickness(:)
92  endif
93 
94  if (present(interp_cs)) cs%interp_CS = interp_cs
95 end subroutine set_hycom_params
96 
97 !> Build a HyCOM coordinate column
98 subroutine build_hycom1_column(CS, eqn_of_state, nz, depth, h, T, S, p_col, &
99  z_col, z_col_new, zScale, h_neglect, h_neglect_edge)
100  type(hycom_cs), intent(in) :: cs !< Coordinate control structure
101  type(eos_type), pointer :: eqn_of_state !< Equation of state structure
102  integer, intent(in) :: nz !< Number of levels
103  real, intent(in) :: depth !< Depth of ocean bottom (positive [H ~> m or kg m-2])
104  real, dimension(nz), intent(in) :: t !< Temperature of column [degC]
105  real, dimension(nz), intent(in) :: s !< Salinity of column [ppt]
106  real, dimension(nz), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
107  real, dimension(nz), intent(in) :: p_col !< Layer pressure [R L2 T-2 ~> Pa]
108  real, dimension(nz+1), intent(in) :: z_col !< Interface positions relative to the surface [H ~> m or kg m-2]
109  real, dimension(CS%nk+1), intent(inout) :: z_col_new !< Absolute positions of interfaces [H ~> m or kg m-2]
110  real, optional, intent(in) :: zscale !< Scaling factor from the input coordinate thicknesses in [Z ~> m]
111  !! to desired units for zInterface, perhaps GV%Z_to_H.
112  real, optional, intent(in) :: h_neglect !< A negligibly small width for the purpose of
113  !! cell reconstruction [H ~> m or kg m-2]
114  real, optional, intent(in) :: h_neglect_edge !< A negligibly small width for the purpose of
115  !! edge value calculation [H ~> m or kg m-2]
116 
117  ! Local variables
118  integer :: k
119  real, dimension(nz) :: rho_col ! Layer densities in a column [R ~> kg m-3]
120  real, dimension(CS%nk) :: h_col_new ! New layer thicknesses
121  real :: z_scale ! A scaling factor from the input thicknesses to the target thicknesses,
122  ! perhaps 1 or a factor in [H Z-1 ~> 1 or kg m-3]
123  real :: stretching ! z* stretching, converts z* to z [nondim].
124  real :: nominal_z ! Nominal depth of interface when using z* [H ~> m or kg m-2]
125  logical :: maximum_depths_set ! If true, the maximum depths of interface have been set.
126  logical :: maximum_h_set ! If true, the maximum layer thicknesses have been set.
127 
128  maximum_depths_set = allocated(cs%max_interface_depths)
129  maximum_h_set = allocated(cs%max_layer_thickness)
130 
131  z_scale = 1.0 ; if (present(zscale)) z_scale = zscale
132 
133  ! Work bottom recording potential density
134  call calculate_density(t, s, p_col, rho_col, eqn_of_state)
135  ! This ensures the potential density profile is monotonic
136  ! although not necessarily single valued.
137  do k = nz-1, 1, -1
138  rho_col(k) = min( rho_col(k), rho_col(k+1) )
139  enddo
140 
141  ! Interpolates for the target interface position with the rho_col profile
142  ! Based on global density profile, interpolate to generate a new grid
143  call build_and_interpolate_grid(cs%interp_CS, rho_col, nz, h(:), z_col, &
144  cs%target_density, cs%nk, h_col_new, z_col_new, h_neglect, h_neglect_edge)
145 
146  ! Sweep down the interfaces and make sure that the interface is at least
147  ! as deep as a nominal target z* grid
148  nominal_z = 0.
149  stretching = z_col(nz+1) / depth ! Stretches z* to z
150  do k = 2, cs%nk+1
151  nominal_z = nominal_z + (z_scale * cs%coordinateResolution(k-1)) * stretching
152  z_col_new(k) = max( z_col_new(k), nominal_z )
153  z_col_new(k) = min( z_col_new(k), z_col(nz+1) )
154  enddo
155 
156  if (maximum_depths_set .and. maximum_h_set) then ; do k=2,cs%nk
157  ! The loop bounds are 2 & nz so the top and bottom interfaces do not move.
158  ! Recall that z_col_new is positive downward.
159  z_col_new(k) = min(z_col_new(k), cs%max_interface_depths(k), &
160  z_col_new(k-1) + cs%max_layer_thickness(k-1))
161  enddo ; elseif (maximum_depths_set) then ; do k=2,cs%nk
162  z_col_new(k) = min(z_col_new(k), cs%max_interface_depths(k))
163  enddo ; elseif (maximum_h_set) then ; do k=2,cs%nk
164  z_col_new(k) = min(z_col_new(k), z_col_new(k-1) + cs%max_layer_thickness(k-1))
165  enddo ; endif
166 end subroutine build_hycom1_column
167 
168 end module coord_hycom
Calculates density of sea water from T, S and P.
Definition: MOM_EOS.F90:68
Regrid columns for the HyCOM coordinate.
Definition: coord_hycom.F90:2
Routines for error handling and I/O management.
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
Control structure for regrid_interp module.
Control structure containing required parameters for the HyCOM coordinate.
Definition: coord_hycom.F90:13
Vertical interpolation for regridding.
A control structure for the equation of state.
Definition: MOM_EOS.F90:108