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)
 
subroutine, public set_hycom_params (CS, max_interface_depths, max_layer_thickness, interp_CS)
 
subroutine, public build_hycom1_column (CS, eqn_of_state, nz, depth, h, T, S, p_col, z_col, z_col_new)
 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(nz+1), intent(inout)  z_col_new 
)

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 in H)
[in]sT and S for column
[in]hLayer thicknesses, in m
[in]p_colLayer pressure in Pa
[in,out]z_col_newAbsolute positions of interfaces

Definition at line 95 of file coord_hycom.F90.

References regrid_interp::build_and_interpolate_grid().

Referenced by mom_regridding::build_grid_hycom1().

95  type(hycom_cs), intent(in) :: cs !< Coordinate control structure
96  type(eos_type), pointer :: eqn_of_state !< Equation of state structure
97  integer, intent(in) :: nz !< Number of levels
98  real, intent(in) :: depth !< Depth of ocean bottom (positive in H)
99  real, dimension(nz), intent(in) :: t, s !< T and S for column
100  real, dimension(nz), intent(in) :: h !< Layer thicknesses, in m
101  real, dimension(nz), intent(in) :: p_col !< Layer pressure in Pa
102  real, dimension(nz+1), intent(in) :: z_col ! Interface positions relative to the surface in H units (m or kg m-2)
103  real, dimension(nz+1), intent(inout) :: z_col_new !< Absolute positions of interfaces
104 
105  ! Local variables
106  integer :: k
107  real, dimension(nz) :: rho_col, h_col_new ! Layer quantities
108  real :: stretching ! z* stretching, converts z* to z.
109  real :: nominal_z ! Nominal depth of interface is using z* (m or Pa)
110  real :: hnew
111  logical :: maximum_depths_set ! If true, the maximum depths of interface have been set.
112  logical :: maximum_h_set ! If true, the maximum layer thicknesses have been set.
113 
114  maximum_depths_set = allocated(cs%max_interface_depths)
115  maximum_h_set = allocated(cs%max_layer_thickness)
116 
117  ! Work bottom recording potential density
118  call calculate_density(t, s, p_col, rho_col, 1, nz, eqn_of_state)
119  ! This ensures the potential density profile is monotonic
120  ! although not necessarily single valued.
121  do k = nz-1, 1, -1
122  rho_col(k) = min( rho_col(k), rho_col(k+1) )
123  enddo
124 
125  ! Interpolates for the target interface position with the rho_col profile
126  ! Based on global density profile, interpolate to generate a new grid
127  call build_and_interpolate_grid(cs%interp_CS, rho_col, nz, h(:), z_col, &
128  cs%target_density, nz, h_col_new, z_col_new)
129 
130  ! Sweep down the interfaces and make sure that the interface is at least
131  ! as deep as a nominal target z* grid
132  nominal_z = 0.
133  stretching = z_col(nz+1) / depth ! Stretches z* to z
134  do k = 2, nz+1
135  nominal_z = nominal_z + cs%coordinateResolution(k-1) * stretching
136  z_col_new(k) = max( z_col_new(k), nominal_z )
137  z_col_new(k) = min( z_col_new(k), z_col(nz+1) )
138  enddo
139 
140  if (maximum_depths_set .and. maximum_h_set) then ; do k=2,nz
141  ! The loop bounds are 2 & nz so the top and bottom interfaces do not move.
142  ! Recall that z_col_new is positive downward.
143  z_col_new(k) = min(z_col_new(k), cs%max_interface_depths(k), &
144  z_col_new(k-1) + cs%max_layer_thickness(k-1))
145  enddo ; elseif (maximum_depths_set) then ; do k=2,nz
146  z_col_new(k) = min(z_col_new(k), cs%max_interface_depths(k))
147  enddo ; elseif (maximum_h_set) then ; do k=2,nz
148  z_col_new(k) = min(z_col_new(k), z_col_new(k-1) + cs%max_layer_thickness(k-1))
149  enddo ; endif
Here is the call graph for this function:
Here is the caller graph for this function:

◆ end_coord_hycom()

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

Definition at line 57 of file coord_hycom.F90.

Referenced by mom_regridding::end_regridding().

57  type(hycom_cs), pointer :: cs
58 
59  ! nothing to do
60  if (.not. associated(cs)) return
61  deallocate(cs%coordinateResolution)
62  deallocate(cs%target_density)
63  if (allocated(cs%max_interface_depths)) deallocate(cs%max_interface_depths)
64  if (allocated(cs%max_layer_thickness)) deallocate(cs%max_layer_thickness)
65  deallocate(cs)
Here is the caller graph for this function:

◆ init_coord_hycom()

subroutine, public coord_hycom::init_coord_hycom ( type(hycom_cs), pointer  CS,
integer, intent(in)  nk,
real, dimension(:), intent(in)  coordinateResolution,
real, dimension(:), 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]coordinateresolutionZ-space thicknesses (m)
[in]target_densityInterface target densities (kg/m3)
[in]interp_csControls for interpolation

Definition at line 39 of file coord_hycom.F90.

References mom_error_handler::mom_error().

39  type(hycom_cs), pointer :: cs !< Unassociated pointer to hold the control structure
40  integer, intent(in) :: nk !< Number of layers in generated grid
41  real, dimension(:), intent(in) :: coordinateresolution !< Z-space thicknesses (m)
42  real, dimension(:), intent(in) :: target_density !< Interface target densities (kg/m3)
43  type(interp_cs_type), intent(in) :: interp_cs !< Controls for interpolation
44 
45  if (associated(cs)) call mom_error(fatal, "init_coord_hycom: CS already associated!")
46  allocate(cs)
47  allocate(cs%coordinateResolution(nk))
48  allocate(cs%target_density(nk+1))
49 
50  cs%nk = nk
51  cs%coordinateResolution(:) = coordinateresolution(:)
52  cs%target_density(:) = target_density(:)
53  cs%interp_CS = interp_cs
Here is the call graph for this function:

◆ 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 
)

Definition at line 69 of file coord_hycom.F90.

References mom_error_handler::mom_error().

69  type(hycom_cs), pointer :: cs
70  real, optional, dimension(:), intent(in) :: max_interface_depths
71  real, optional, dimension(:), intent(in) :: max_layer_thickness
72  type(interp_cs_type), optional, intent(in) :: interp_cs
73 
74  if (.not. associated(cs)) call mom_error(fatal, "set_hycom_params: CS not associated")
75 
76  if (present(max_interface_depths)) then
77  if (size(max_interface_depths) /= cs%nk+1) &
78  call mom_error(fatal, "set_hycom_params: max_interface_depths inconsistent size")
79  allocate(cs%max_interface_depths(cs%nk+1))
80  cs%max_interface_depths(:) = max_interface_depths(:)
81  endif
82 
83  if (present(max_layer_thickness)) then
84  if (size(max_layer_thickness) /= cs%nk) &
85  call mom_error(fatal, "set_hycom_params: max_layer_thickness inconsistent size")
86  allocate(cs%max_layer_thickness(cs%nk))
87  cs%max_layer_thickness(:) = max_layer_thickness(:)
88  endif
89 
90  if (present(interp_cs)) cs%interp_CS = interp_cs
Here is the call graph for this function: