MOM6
coord_hycom Module Reference

Detailed Description

Regrid columns for the HyCOM coordinate.

Data Types

type  hycom_cs
 Control structure containing required parameters for the HyCOM coordinate. More...
 

Functions/Subroutines

subroutine, public init_coord_hycom (CS, nk, coordinateResolution, target_density, interp_CS)
 Initialise a hycom_CS with pointers to parameters. More...
 
subroutine, public end_coord_hycom (CS)
 This subroutine deallocates memory in the control structure for the coord_hycom module. More...
 
subroutine, public set_hycom_params (CS, max_interface_depths, max_layer_thickness, interp_CS)
 This subroutine can be used to set the parameters for the coord_hycom module. More...
 
subroutine, public build_hycom1_column (CS, eqn_of_state, nz, depth, h, T, S, p_col, z_col, z_col_new, zScale, h_neglect, h_neglect_edge)
 Build a HyCOM coordinate column. More...
 

Function/Subroutine Documentation

◆ build_hycom1_column()

subroutine, public coord_hycom::build_hycom1_column ( type(hycom_cs), intent(in)  CS,
type(eos_type), pointer  eqn_of_state,
integer, intent(in)  nz,
real, intent(in)  depth,
real, dimension(nz), intent(in)  h,
real, dimension(nz), intent(in)  T,
real, dimension(nz), intent(in)  S,
real, dimension(nz), intent(in)  p_col,
real, dimension(nz+1), intent(in)  z_col,
real, dimension(cs%nk+1), intent(inout)  z_col_new,
real, intent(in), optional  zScale,
real, intent(in), optional  h_neglect,
real, intent(in), optional  h_neglect_edge 
)

Build a HyCOM coordinate column.

Parameters
[in]csCoordinate control structure
eqn_of_stateEquation of state structure
[in]nzNumber of levels
[in]depthDepth of ocean bottom (positive [H ~> m or kg m-2])
[in]tTemperature of column [degC]
[in]sSalinity of column [ppt]
[in]hLayer thicknesses, in [m] or [H ~> m or kg m-2]
[in]p_colLayer pressure [Pa]
[in]z_colInterface positions relative to the surface [H ~> m or kg m-2]
[in,out]z_col_newAbsolute positions of interfaces
[in]zscaleScaling factor from the input thicknesses in [m] to desired units for zInterface, perhaps m_to_H.
[in]h_neglectA negligibly small width for the purpose of cell reconstructions in the same units as h.
[in]h_neglect_edgeA negligibly small width for the purpose of edge value calculations in the same units as h0.

Definition at line 99 of file coord_hycom.F90.

99  type(hycom_CS), intent(in) :: CS !< Coordinate control structure
100  type(EOS_type), pointer :: eqn_of_state !< Equation of state structure
101  integer, intent(in) :: nz !< Number of levels
102  real, intent(in) :: depth !< Depth of ocean bottom (positive [H ~> m or kg m-2])
103  real, dimension(nz), intent(in) :: T !< Temperature of column [degC]
104  real, dimension(nz), intent(in) :: S !< Salinity of column [ppt]
105  real, dimension(nz), intent(in) :: h !< Layer thicknesses, in [m] or [H ~> m or kg m-2]
106  real, dimension(nz), intent(in) :: p_col !< Layer pressure [Pa]
107  real, dimension(nz+1), intent(in) :: z_col !< Interface positions relative to the surface [H ~> m or kg m-2]
108  real, dimension(CS%nk+1), intent(inout) :: z_col_new !< Absolute positions of interfaces
109  real, optional, intent(in) :: zScale !< Scaling factor from the input thicknesses in [m]
110  !! to desired units for zInterface, perhaps m_to_H.
111  real, optional, intent(in) :: h_neglect !< A negligibly small width for the
112  !! purpose of cell reconstructions
113  !! in the same units as h.
114  real, optional, intent(in) :: h_neglect_edge !< A negligibly small width
115  !! for the purpose of edge value calculations
116  !! in the same units as h0.
117 
118  ! Local variables
119  integer :: k
120  real, dimension(nz) :: rho_col ! Layer quantities
121  real, dimension(CS%nk) :: h_col_new ! New layer thicknesses
122  real :: z_scale
123  real :: stretching ! z* stretching, converts z* to z.
124  real :: nominal_z ! Nominal depth of interface when using z* [Z ~> m]
125  real :: hNew
126  logical :: maximum_depths_set ! If true, the maximum depths of interface have been set.
127  logical :: maximum_h_set ! If true, the maximum layer thicknesses have been set.
128 
129  maximum_depths_set = allocated(cs%max_interface_depths)
130  maximum_h_set = allocated(cs%max_layer_thickness)
131 
132  z_scale = 1.0 ; if (present(zscale)) z_scale = zscale
133 
134  ! Work bottom recording potential density
135  call calculate_density(t, s, p_col, rho_col, 1, nz, eqn_of_state)
136  ! This ensures the potential density profile is monotonic
137  ! although not necessarily single valued.
138  do k = nz-1, 1, -1
139  rho_col(k) = min( rho_col(k), rho_col(k+1) )
140  enddo
141 
142  ! Interpolates for the target interface position with the rho_col profile
143  ! Based on global density profile, interpolate to generate a new grid
144  call build_and_interpolate_grid(cs%interp_CS, rho_col, nz, h(:), z_col, &
145  cs%target_density, cs%nk, h_col_new, z_col_new, h_neglect, h_neglect_edge)
146 
147  ! Sweep down the interfaces and make sure that the interface is at least
148  ! as deep as a nominal target z* grid
149  nominal_z = 0.
150  stretching = z_col(nz+1) / depth ! Stretches z* to z
151  do k = 2, cs%nk+1
152  nominal_z = nominal_z + (z_scale * cs%coordinateResolution(k-1)) * stretching
153  z_col_new(k) = max( z_col_new(k), nominal_z )
154  z_col_new(k) = min( z_col_new(k), z_col(nz+1) )
155  enddo
156 
157  if (maximum_depths_set .and. maximum_h_set) then ; do k=2,cs%nk
158  ! The loop bounds are 2 & nz so the top and bottom interfaces do not move.
159  ! Recall that z_col_new is positive downward.
160  z_col_new(k) = min(z_col_new(k), cs%max_interface_depths(k), &
161  z_col_new(k-1) + cs%max_layer_thickness(k-1))
162  enddo ; elseif (maximum_depths_set) then ; do k=2,cs%nk
163  z_col_new(k) = min(z_col_new(k), cs%max_interface_depths(k))
164  enddo ; elseif (maximum_h_set) then ; do k=2,cs%nk
165  z_col_new(k) = min(z_col_new(k), z_col_new(k-1) + cs%max_layer_thickness(k-1))
166  enddo ; endif

◆ end_coord_hycom()

subroutine, public coord_hycom::end_coord_hycom ( type(hycom_cs), pointer  CS)

This subroutine deallocates memory in the control structure for the coord_hycom module.

Parameters
csCoordinate control structure

Definition at line 59 of file coord_hycom.F90.

59  type(hycom_CS), pointer :: CS !< Coordinate control structure
60 
61  ! nothing to do
62  if (.not. associated(cs)) return
63  deallocate(cs%coordinateResolution)
64  deallocate(cs%target_density)
65  if (allocated(cs%max_interface_depths)) deallocate(cs%max_interface_depths)
66  if (allocated(cs%max_layer_thickness)) deallocate(cs%max_layer_thickness)
67  deallocate(cs)

◆ init_coord_hycom()

subroutine, public coord_hycom::init_coord_hycom ( type(hycom_cs), pointer  CS,
integer, intent(in)  nk,
real, dimension(nk), intent(in)  coordinateResolution,
real, dimension(nk+1), intent(in)  target_density,
type(interp_cs_type), intent(in)  interp_CS 
)

Initialise a hycom_CS with pointers to parameters.

Parameters
csUnassociated pointer to hold the control structure
[in]nkNumber of layers in generated grid
[in]coordinateresolutionNominal near-surface resolution [m]
[in]target_densityInterface target densities [kg m-3]
[in]interp_csControls for interpolation

Definition at line 40 of file coord_hycom.F90.

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 [m]
43  real, dimension(nk+1),intent(in) :: target_density !< Interface target densities [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

◆ set_hycom_params()

subroutine, public coord_hycom::set_hycom_params ( type(hycom_cs), pointer  CS,
real, dimension(:), intent(in), optional  max_interface_depths,
real, dimension(:), intent(in), optional  max_layer_thickness,
type(interp_cs_type), intent(in), optional  interp_CS 
)

This subroutine can be used to set the parameters for the coord_hycom module.

Parameters
csCoordinate control structure
[in]max_interface_depthsMaximum depths of interfaces in m
[in]max_layer_thicknessMaximum thicknesses of layers in m
[in]interp_csControls for interpolation

Definition at line 72 of file coord_hycom.F90.

72  type(hycom_CS), pointer :: CS !< Coordinate control structure
73  real, dimension(:), optional, intent(in) :: max_interface_depths !< Maximum depths of interfaces in m
74  real, dimension(:), optional, intent(in) :: max_layer_thickness !< Maximum thicknesses of layers in m
75  type(interp_CS_type), optional, intent(in) :: interp_CS !< Controls for interpolation
76 
77  if (.not. associated(cs)) call mom_error(fatal, "set_hycom_params: CS not associated")
78 
79  if (present(max_interface_depths)) then
80  if (size(max_interface_depths) /= cs%nk+1) &
81  call mom_error(fatal, "set_hycom_params: max_interface_depths inconsistent size")
82  allocate(cs%max_interface_depths(cs%nk+1))
83  cs%max_interface_depths(:) = max_interface_depths(:)
84  endif
85 
86  if (present(max_layer_thickness)) then
87  if (size(max_layer_thickness) /= cs%nk) &
88  call mom_error(fatal, "set_hycom_params: max_layer_thickness inconsistent size")
89  allocate(cs%max_layer_thickness(cs%nk))
90  cs%max_layer_thickness(:) = max_layer_thickness(:)
91  endif
92 
93  if (present(interp_cs)) cs%interp_CS = interp_cs